SplineReg {GeDS} | R Documentation |
Estimation of the coefficients of a predictor model with spline and possibly parametric components.
Description
Functions that estimate the coefficients of a predictor model involving a spline component and possibly a parametric component applying (Iteratively Re-weighted) Least Squares (IR)LS iteration.
Usage
SplineReg_LM(
X,
Y,
Z = NULL,
offset = rep(0, length(X)),
weights = rep(1, length(X)),
InterKnots,
n,
extr = range(X),
prob = 0.95,
coefficients = NULL,
only_predictions = FALSE
)
SplineReg_GLM(
X,
Y,
Z,
offset = rep(0, nobs),
weights = rep(1, length(X)),
InterKnots,
n,
extr = range(X),
family,
mustart,
inits = NULL,
etastart = NULL
)
Arguments
X |
a numeric vector containing |
Y |
a vector of size |
Z |
a design matrix with |
offset |
a vector of size |
weights |
an optional vector of ‘prior weights’ to be put on the observations in the fitting process in case the user requires weighted fitting. It is a vector of 1s by default. |
InterKnots |
a numeric vector containing the locations of the internal knots necessary to compute the B-splines. In GeDS these are the internal knots in a current iteration of stage A. |
n |
integer value specifying the order of the spline to be evaluated.
It should be 2 (linear spline), 3 (quadratic spline) or 4 (cubic spline).
Non-integer values will be passed to the function |
extr |
optional numeric vector of 2 elements representing the left-most
and right-most limits of the interval embedding the sample values of
|
prob |
the confidence level to be used for the confidence bands in the
|
coefficients |
optional vector of spline coefficients. If provided,
|
only_predictions |
logical, if |
family |
a description of the error distribution and link function to be
used in the model. This can be a character string naming a family function, a
family function or the result of a call to a family function. See
|
mustart |
initial values for the vector of means in the IRLS estimation.
Must be a vector of length |
inits |
a numeric vector of length
|
etastart |
initial values for the predictor in the IRLS estimation
(alternative to providing either |
Details
The functions estimate the coefficients of a predictor model with a spline
component (and possibly a parametric component) for a given, fixed order and
vector of knots of the spline and a specified distribution of the response
variable (from the Exponential Family). The functions SplineReg_LM
and
SplineReg_GLM
are based correspondingly on LS and IRLS and used
correspondingly in NGeDS
and GGeDS
, to estimate
the coefficients of the final GeDS fits of stage B, after their knots have
been positioned to coincide with the Greville abscissas of the knots of the
linear fit from stage A (see Dimitrova et al. 2023). Additional inference
related quantities are also computed (see Value below). The function
SplineReg_GLM
is also used to estimate the coefficients of the linear
GeDS fit of stage A within GGeDS
, whereas in
NGeDS
this estimation is performed internally leading to faster
R code.
In addition SplineReg_LM
computes some useful quantities among which
confidence intervals and the Control Polygon (see Section 2 of
Kaishev et al. 2016).
The confidence intervals contained in the output slot NCI
are
approximate local bands obtained considering the knots as fixed constants.
Hence the columns of the design matrix are seen as covariates and standard
methodology relying on the se.fit
option of predict.lm
or
predict.glm
is used. In the ACI
slot, asymptotic confidence
intervals are provided, following Kaishev et al (2006).
As mentioned, SplineReg_GLM
is intensively used in Stage A of the GeDS
algorithm implemented in GGeDS
and in order to make it as fast
as possible input data validation is mild. Hence it is expected that the user
checks carefully the input parameters before using SplineReg_GLM
. The
"Residuals
" in the output of this function are similar to the so
called “working residuals" in the glm
function.
"Residuals
" are the residuals r_i
used in the knot placement
procedure, i.e.
r_i= (y_i - \hat{\mu}_i){d \mu_i \over d \eta_i },
but
in contrast to glm
“working residuals", they are
computed using the final IRLS fitted \hat{\mu}_i
. "Residuals
"
are then used in locating the knots of the linear spline fit of Stage A.
In SplineReg_GLM
confidence intervals are not computed.
Value
A list
containing:
Theta |
a vector containing the fitted coefficients of the spline regression component and the parametric component of the predictor model. |
Predicted |
a vector of |
Residuals |
a vector containing the normal regression residuals if
|
RSS |
the deviance for the fitted predictor model, defined as in
Dimitrova et al. (2023), which for |
NCI |
a list containing the lower ( |
Basis |
the matrix of B-spline regression functions and the covariates of the parametric part evaluated at the sample values of the covariate(s). |
Polygon |
a list containing x-y coordinates (" |
deviance |
a vector containing deviances computed at each IRLS step
(computed only with the |
temporary |
the result of the function |
References
Kaishev, V. K., Dimitrova, D. S., Haberman, S. & Verrall, R. J. (2006).
Geometrically designed, variable know regression splines: asymptotics and
inference (Statistical Research Paper No. 28).
London, UK: Faculty of Actuarial Science & Insurance, City University London.
URL: openaccess.city.ac.uk
Kaishev, V.K., Dimitrova, D.S., Haberman, S., & Verrall, R.J. (2016).
Geometrically designed, variable knot regression splines.
Computational Statistics, 31, 1079–1105.
DOI: doi:10.1007/s00180-015-0621-7
Dimitrova, D. S., Kaishev, V. K., Lattuada, A. and Verrall, R. J. (2023).
Geometrically designed variable knot splines in generalized (non-)linear
models. Applied Mathematics and Computation, 436.
DOI: doi:10.1016/j.amc.2022.127493
See Also
NGeDS
, GGeDS
, Fitters
,
IRLSfit
, lm
and
glm.fit
.