MGCG_smooth: High-dimensional spline smoothing using a matrix-free...

Description Usage Arguments Value References Examples

View source: R/mgss_main.R

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

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
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 nrow(X) as the variable of interest.

w

Damping factor of the Jacobi smoother. Defaults to 0.1.

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 prod(m+q+1) as starting value for the MGCG-method. Defaults to zero.

K_max

Positive integer as upper bound for the number of MGCG-iterations. Defaults to prod(m+q+1).

tolerance

Positive number as error tolerance for the stopping criterion of the MGCG-method. Defaults to 1e-6.

print_error

Logical, indicating if the iteration error should be printed or not.

coarse_grid_solver

Utilized coarse grid solver. Either "PCG" for diagonal preconditioned CG or "Cholesky" for Cholesky decomposition. Defaults to "Cholesky".

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

1
2
3
4
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)

mgss documentation built on May 10, 2021, 9:07 a.m.