MGCG_smooth {mgss} | R Documentation |
High-dimensional spline smoothing using a matrix-free multigrid preconditioned CG-method.
Description
Fits a smooth spline to a set of given observations using penalized splines with curvature penalty and multiple covariates. The underlying linear system is solved with a matrix-free preconditioned conjugated gradient method using a geometric multigrid method as preconditioner.
Usage
MGCG_smooth(
G,
q,
lambda,
X,
y,
w = 0.1,
nu = c(3, 1),
alpha_start = NULL,
K_max = NULL,
tolerance = 1e-06,
print_error = TRUE,
coarse_grid_solver = "Cholesky"
)
Arguments
G |
Positive integer greater than one for the maximum number of grids. |
q |
Vector of positive integers. Each entry gives the spline degree for the respective covariate. |
lambda |
Positive number as weight for the penalty term. |
X |
Matrix containing the covariates as columns and the units as rows. |
y |
Vector of length |
w |
Damping factor of the Jacobi smoother. Defaults to |
nu |
Two-dimensional vector of non-negative integers. Gives the number of pre- and post-smoothing steps in the multigrid algorithm. |
alpha_start |
Vector of length |
K_max |
Positive integer as upper bound for the number of MGCG-iterations. Defaults to |
tolerance |
Positive number as error tolerance for the stopping criterion of the MGCG-method. Defaults to |
print_error |
Logical, indicating if the iteration error should be printed or not. |
coarse_grid_solver |
Utilized coarse grid solver. Either |
Value
Returns a list containing the input m = 2^G-1
, q
, and Omega
. Further gives the fitted spline coefficients alpha
, the fitted values fitted_values
, the residuals residuals
, the root mean squared error rmse
and the R-squared value R_squared
.
References
Siebenborn, M. and Wagner, J. (2019) A Multigrid Preconditioner for Tensor Product Spline Smoothing. arXiv:1901.00654
Examples
data <- generate_test_data(100, 2)
X <- data$X_train
y <- data$y_train
MGCG_smooth(G = 3, q = c(3,3), lambda = 0.1, w = 0.8, X = X, y = y)