slgp | R Documentation |
This function builds and trains an SLGP model based on a specified formula and data. The SLGP is a finite-rank Gaussian process model for conditional density estimation, trained using MAP, MCMC, Laplace approximation, or left untrained ("none").
slgp(
formula,
data,
epsilonStart = NULL,
method,
basisFunctionsUsed,
interpolateBasisFun = "NN",
nIntegral = 51,
nDiscret = 51,
hyperparams = NULL,
predictorsUpper = NULL,
predictorsLower = NULL,
responseRange = NULL,
sigmaEstimationMethod = "none",
seed = NULL,
opts_BasisFun = list(),
BasisFunParam = NULL,
opts = list(),
verbose = FALSE
)
formula |
A formula specifying the model structure, with the response on the left-hand side and covariates on the right. |
data |
A data frame containing the variables used in the formula. |
epsilonStart |
Optional numeric vector of initial weights for the finite-rank GP:
|
method |
Character string specifying the training method: one of {"none", "MCMC", "MAP", "Laplace"}. |
basisFunctionsUsed |
Character string describing the basis function type: one of "inducing points", "RFF", "Discrete FF", "filling FF", or "custom cosines". |
interpolateBasisFun |
Character string indicating how to evaluate basis functions: "nothing" (exact eval), "NN" (nearest-neighbor), or "WNN" (weighted inverse-distance). Default is "NN". |
nIntegral |
Number of quadrature points used for numerical integration over the response domain. |
nDiscret |
Integer controlling the resolution of the interpolation grid (used only for "NN" or "WNN"). |
hyperparams |
Optional list of hyperparameters. Should contain:
|
predictorsUpper |
Optional numeric vector for the upper bounds of the covariates (used for scaling). |
predictorsLower |
Optional numeric vector for the lower bounds of the covariates. |
responseRange |
Optional numeric vector of length 2 with the lower and upper bounds of the response. |
sigmaEstimationMethod |
Method to heuristically estimate the variance |
seed |
Optional integer for reproducibility. |
opts_BasisFun |
List of optional configuration parameters passed to the basis function initializer. |
BasisFunParam |
Optional list of precomputed basis function parameters. |
opts |
Optional list of extra settings passed to inference routines (e.g., |
verbose |
Logical; if |
An object of S4 class SLGP-class
, containing:
Matrix of posterior (or prior) draws of the SLGP coefficients \epsilon_i
.
List of fitted or provided hyperparameters.
Log-posterior (if MAP or Laplace used).
Estimation method used.
Other internal information such as ranges, basis settings, and data.
Gautier, Athénaïs (2023). "Modelling and Predicting Distribution-Valued Fields with Applications to Inversion Under Uncertainty." Thesis, Universität Bern, Bern. https://boristheses.unibe.ch/4377/
# Load Boston housing dataset
library(MASS)
data("Boston")
# Set input and output ranges manually (you can also use range(Boston$age), etc.)
range_x <- c(0, 100)
range_response <- c(0, 50)
#' #Create a SLGP model but don't fit it
modelPrior <- slgp(medv ~ age, # Use a formula to specify response and covariates
data = Boston, # Use the original Boston housing data
method = "none", # No training
basisFunctionsUsed = "RFF", # Random Fourier Features
sigmaEstimationMethod = "heuristic", # Auto-tune sigma2 (more stable)
predictorsLower = range_x[1], # Lower bound for 'age'
predictorsUpper = range_x[2], # Upper bound for 'age'
responseRange = range_response, # Range for 'medv'
opts_BasisFun = list(nFreq = 200, # Use 200 Fourier features
MatParam = 5/2), # Matern 5/2 kernel
seed = 1) # Reproducibility
# Train an SLGP model using MAP estimation and RFF basis
modelMAP <- slgp(medv ~ age, # Use a formula to specify response and covariates
data = Boston, # Use the original Boston housing data
method = "MAP", # Train using Maximum A Posteriori estimation
basisFunctionsUsed = "RFF", # Random Fourier Features
sigmaEstimationMethod = "heuristic", # Auto-tune sigma2 (more stable)
predictorsLower = range_x[1], # Lower bound for 'age'
predictorsUpper = range_x[2], # Upper bound for 'age'
responseRange = range_response, # Range for 'medv'
opts_BasisFun = list(nFreq = 200, # Use 200 Fourier features
MatParam = 5/2), # Matern 5/2 kernel
seed = 1) # Reproducibility
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.