File: naturalSpline.R

package info (click to toggle)
r-cran-splines2 0.4.1-1
  • links: PTS, VCS
  • area: main
  • in suites: bullseye, sid
  • size: 1,512 kB
  • sloc: cpp: 1,988; ansic: 165; sh: 13; makefile: 2
file content (123 lines) | stat: -rw-r--r-- 4,695 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
##
## R package splines2 by Wenjie Wang and Jun Yan
## Copyright (C) 2016-2021
##
## This file is part of the R package splines2.
##
## The R package splines2 is free software: You can redistribute it and/or
## modify it under the terms of the GNU General Public License as published by
## the Free Software Foundation, either version 3 of the License, or any later
## version (at your option). See the GNU General Public License at
## <https://www.gnu.org/licenses/> for details.
##
## The R package splines2 is distributed in the hope that it will be useful,
## but WITHOUT ANY WARRANTY without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
##

##' Natural Cubic Spline Basis for Polynomial Splines
##'
##' Generates the nonnegative natural cubic spline basis matrix, the
##' corresponding integrals (from the left boundary knot), or derivatives of
##' given order.  Each basis is assumed to follow a linear trend for \code{x}
##' outside of boundary.
##'
##' It is an implementation of the natural spline basis based on B-spline basis,
##' which utilizes the close-form null space that can be derived from the
##' recursive formula for the second derivatives of B-splines.  The constructed
##' spline bases are intended to be nonnegative within boundary with second
##' derivatives being zeros at boundary knots.
##'
##' A similar implementation is provided by \code{splines::ns}, which uses QR
##' decomposition to find the null space of the second derivatives of B-spline
##' basis at boundary knots.  However, there is no guarantee that the resulting
##' bases are nonnegative within boundary.
##'
##' @inheritParams bSpline
##'
##' @param df Degree of freedom that equals to the column number of returned
##'     matrix.  One can specify \code{df} rather than \code{knots}, then the
##'     function chooses \code{df - 1 - as.integer(intercept)} internal knots at
##'     suitable quantiles of \code{x} ignoring missing values and those
##'     \code{x} outside of the boundary.  Thus, \code{df} must be greater than
##'     or equal to \code{2}.  If internal knots are specified via \code{knots},
##'     the specified \code{df} will be ignored.
##' @param derivs A nonnegative integer specifying the order of derivatives of
##'     natural splines. The default value is \code{0L} for the spline bases.
##' @param integral A logical value.  The default value is \code{FALSE}.  If
##'     \code{TRUE}, this function will return the integrated natural splines
##'     from the left boundary knot.
##'
##' @return A numeric matrix of \code{length(x)} rows and \code{df}
##'     columns if \code{df} is specified or \code{length(knots) + 1 +
##'     as.integer(intercept)} columns if \code{knots} are specified instead.
##'     Attributes that correspond to the arguments specified are returned for
##'     usage of other functions in this package.
##'
##' @example inst/examples/ex-naturalSpline.R
##'
##' @seealso
##' \code{\link{bSpline}} for B-splines;
##' \code{\link{mSpline}} for M-splines;
##' \code{\link{iSpline}} for I-splines.
##'
##' @export
naturalSpline <- function(x, df = NULL, knots = NULL,
                          intercept = FALSE, Boundary.knots = NULL,
                          derivs = 0L, integral = FALSE,
                          ...)
{
    ## check inputs
    if ((derivs <- as.integer(derivs)) < 0) {
        stop("The 'derivs' must be a nonnegative integer.")
    }
    if (is.null(df)) {
        df <- 0L
    } else {
        df <- as.integer(df)
        if (df < 2) {
            stop("The 'df' must be >= 2.")
        }
    }
    knots <- null2num0(knots)
    Boundary.knots <- null2num0(Boundary.knots)
    ## take care of possible NA's in `x`
    nax <- is.na(x)
    if (all(nax)) {
        stop("The 'x' cannot be all NA's!")
    }
    ## remove NA's in x
    xx <- if (nas <- any(nax)) {
              x[! nax]
          } else {
              x
          }
    ## call the engine function
    out <- rcpp_naturalSpline(
        x = xx,
        df = df,
        internal_knots = knots,
        boundary_knots = Boundary.knots,
        derivs = derivs,
        integral = integral,
        complete_basis = intercept
    )
    ## keep NA's as is
    if (nas) {
        nmat <- matrix(NA, length(nax), ncol(out))
        nmat[! nax, ] <- out
        saved_attr <- attributes(out)
        saved_attr$dim[1] <- length(nax)
        out <- nmat
        attributes(out) <- saved_attr
        attr(out, "x") <- x
    }
    ## add dimnames for consistency
    name_x <- names(x)
    if (! is.null(name_x)) {
        row.names(out) <- name_x
    }
    ## add class
    class(out) <- c("matrix", "naturalSpline")
    out
}