API
Splinter.BSplineBasisSplinter.ISplineBasisSplinter.MSplineBasisSplinter.NSplineBasisSplinter.bsSplinter.isSplinter.msSplinter.nsSplinter.spline_args
Splinter.BSplineBasis — Method
BSplineBasis(x::AbstractVector{T};
boundary_knots::Union{Tuple{T,T},Nothing}=nothing,
interior_knots::Union{AbstractVector{T},Nothing}=nothing,
order::Int=4,
intercept::Bool=false,
df::Int=order - 1 + Int(intercept),
knots::Union{AbstractVector{T},Nothing}=nothing) where {T<:Real}Calculate a basis for B-splines and return a callable type for evaluating points in that basis.
Keyword Arguments
boundary_knots: boundary knotsinterior_knots: interior knotsorder: order of the splineintercept: bool for whether to include an interceptdf: degrees of freedomknots: full set of knots (excluding repeats)
Keyword Arguments of Returned Callable
center: value to center the splinesderivs: derivatives of the splines
Examples
julia> spline = BSplineBasis(collect(0.0:0.2:1.0); df=3);
julia> spline(collect(0.0:0.2:1.0))
6×3 Array{Float64,2}:
0.0 0.0 0.0
0.384 0.096 0.008
0.432 0.288 0.064
0.288 0.432 0.216
0.096 0.384 0.512
0.0 0.0 1.0Splinter.ISplineBasis — Method
ISplineBasis(x::AbstractVector{T};
boundary_knots::Union{Tuple{T,T},Nothing}=nothing,
interior_knots::Union{AbstractVector{T},Nothing}=nothing,
order::Int=4,
intercept::Bool=false,
df::Int=order + Int(intercept),
knots::Union{AbstractVector{T},Nothing}=nothing) where {T<:Real}Calculate a basis for I-splines and return a callable type for evaluating points in that basis.
Keyword Arguments
boundary_knots: boundary knotsinterior_knots: interior knotsorder: order of the splineintercept: bool for whether to include an intercept (column of ones). This behaviour is different to the splines2 package from R, where intercept=FALSE will drop the first spline term.df: degrees of freedomknots: full set of knots (excluding repeats)
Keyword Arguments of Returned Callable
derivs: derivatives of the splines
Examples
julia> spline = ISplineBasis(collect(0.0:0.2:1.0); df=3);
julia> spline(collect(0.0:0.2:1.0))
6×3 Array{Float64,2}:
0.0 0.0 0.0
0.1808 0.0272 0.0016
0.5248 0.1792 0.0256
0.8208 0.4752 0.1296
0.9728 0.8192 0.4096
1.0 1.0 1.0Splinter.MSplineBasis — Method
MSplineBasis(x::AbstractVector{T};
boundary_knots::Union{Tuple{T,T},Nothing}=nothing,
interior_knots::Union{AbstractVector{T},Nothing}=nothing,
order::Int=4,
intercept::Bool=false,
df::Int=order - 1 + Int(intercept),
knots::Union{AbstractVector{T},Nothing}=nothing) where {T<:Real}Calculate a basis for M-splines and return a callable type for evaluating points in that basis.
Keyword Arguments
boundary_knots: boundary knotsinterior_knots: interior knotsorder: order of the splineintercept: bool for whether to include an intercept (column of ones). This behaviour is different to the Splinter package from R, where intercept=FALSE will drop the first spline term.df: degrees of freedomknots: full set of knots (excluding repeats)
Keyword Arguments of Returned Callable
center: value to center the splinesderivs: derivatives of the splines
Examples
julia> spline = MSplineBasis(collect(0.0:0.2:1.0); df=3);
julia> spline(collect(0.0:0.2:1.0))
6×3 Array{Float64,2}:
0.0 0.0 0.0
1.536 0.384 0.032
1.728 1.152 0.256
1.152 1.728 0.864
0.384 1.536 2.048
0.0 0.0 4.0Splinter.NSplineBasis — Method
NSplineBasis(x::AbstractVector{T};
boundary_knots::Union{Tuple{T,T},Nothing}=nothing,
interior_knots::Union{AbstractVector{T},Nothing}=nothing,
order::Int=4,
intercept::Bool=false,
df::Int=order - 3 + Int(intercept),
knots::Union{AbstractVector{T},Nothing}=nothing) where {T<:Real}Calculate a basis for natural B-splines and return a callable type for evaluating points in that basis.
Keyword Arguments
boundary_knots: boundary knotsinterior_knots: interior knotsorder: order of the splineintercept: bool for whether to include an interceptdf: degrees of freedomknots: full set of knots (excluding repeats)
Keyword Arguments of Returned Callable
center: value to center the splinesderivs: derivatives of the splines
Examples
julia> spline = NSplineBasis(x; df=3);
julia> spline(0.0:0.2:1.0)
6×3 Matrix{Float64}:
0.0 0.0 0.0
-0.100444 0.409332 -0.272888
0.102383 0.540852 -0.359235
0.501759 0.386722 -0.172481
0.418872 0.327383 0.217745
-0.142857 0.428571 0.714286Splinter.bs — Method
bs(x::AbstractVector{T}; center::Union{T,Nothing}=nothing, derivs::Int=0, kwargs...) where {T <: Real}
Calculate a basis for B-splines and return x expressed in that basis.
Keyword arguments
center: value to center the splinesderivs: derivatives of the splines- Further keyword arguments are passed to
BSplineBasis.
Arguments
boundary_knots :: Union{Tuple{T,T},Nothing} = nothing: boundary knotsinterior_knots :: Union{Array{T,1},Nothing} = nothing: interior knotsorder :: Int = 4: order of the splineintercept :: Bool = false: bool for whether to include an interceptdf :: Int = order - 1 + Int(intercept): degrees of freedomknots :: Union{Array{T,1}, Nothing} = nothing: full set of knots (excluding repeats)
Examples
julia> bs(collect(0.0:0.2:1.0), df=3)
6×3 Array{Float64,2}:
0.0 0.0 0.0
0.384 0.096 0.008
0.432 0.288 0.064
0.288 0.432 0.216
0.096 0.384 0.512
0.0 0.0 1.0Splinter.is — Method
is(x::AbstractVector{T}; derivs::Int=0, kwargs...) where {T <: Real}
Calculate a basis for I-splines and return x expressed in that basis.
Keyword arguments
derivs: derivatives of the splines- Further keyword arguments are passed to
ISplineBasis.
Arguments
boundary_knots :: Union{Tuple{T,T},Nothing} = nothing: boundary knotsinterior_knots :: Union{Array{T,1},Nothing} = nothing: interior knotsorder :: Int = 4: order of the splineintercept :: Bool = false: bool for whether to include an interceptdf :: Int = order + Int(intercept): degrees of freedomknots :: Union{Array{T,1}, Nothing} = nothing: full set of knots (excluding repeats)
Examples
julia> is(collect(0.0:0.2:1.0), df=3)
6×3 Array{Float64,2}:
0.0 0.0 0.0
0.1808 0.0272 0.0016
0.5248 0.1792 0.0256
0.8208 0.4752 0.1296
0.9728 0.8192 0.4096
1.0 1.0 1.0Splinter.ms — Method
ms(x::AbstractVector{T}; center::Union{T,Nothing}=nothing, derivs::Int=0, kwargs...) where {T <: Real}
Calculate a basis for M-splines and return x expressed in that basis.
Keyword arguments
center: value to center the splinesderivs: derivatives of the splines- Further keyword arguments are passed to
MSplineBasis.
Arguments
boundary_knots :: Union{Tuple{T,T},Nothing} = nothing: boundary knotsinterior_knots :: Union{Array{T,1},Nothing} = nothing: interior knotsorder :: Int = 4: order of the splineintercept :: Bool = false: bool for whether to include an interceptdf :: Int = order - 1 + Int(intercept): degrees of freedomknots :: Union{Array{T,1}, Nothing} = nothing: full set of knots (excluding repeats)
Examples
julia> ms(collect(0.0:0.2:1.0), df=3)
6×3 Array{Float64,2}:
0.0 0.0 0.0
1.536 0.384 0.032
1.728 1.152 0.256
1.152 1.728 0.864
0.384 1.536 2.048
0.0 0.0 4.0Splinter.ns — Method
ns(x::AbstractVector{T}; center::Union{T,Nothing}=nothing, derivs::Int=0, kwargs...) where {T <: Real}
Calculate a basis for natural B-splines and return x expressed in that basis.
Keyword arguments
center: value to center the splinesderivs: derivatives of the splines- Further keyword arguments are passed to
NSplineBasis.
Arguments
boundary_knots :: Union{Tuple{T,T},Nothing} = nothing: boundary knotsinterior_knots :: Union{Array{T,1},Nothing} = nothing: interior knotsorder :: Int = 4: order of the splineintercept :: Bool = false: bool for whether to include an interceptdf :: Int = order - 1 + Int(intercept): degrees of freedomknots :: Union{Array{T,1}, Nothing} = nothing: full set of knots (excluding repeats)
Examples
julia> ns(0.0:0.2:1.0; df=3)
6×3 Matrix{Float64}:
0.0 0.0 0.0
-0.100444 0.409332 -0.272888
0.102383 0.540852 -0.359235
0.501759 0.386722 -0.172481
0.418872 0.327383 0.217745
-0.142857 0.428571 0.714286Splinter.spline_args — Method
spline_args(x::AbstractVector{T};
boundary_knots::Union{Tuple{T,T},Nothing}=nothing,
interior_knots::Union{Vector{T},Nothing}=nothing,
order::Int=4,
intercept::Bool=false,
df::Int=3 + Int(intercept),
knots::Union{Vector{T},Nothing}=nothing,
knots_offset::Int=0) where {T<:Real}Utility function for processing the spline arguments.
Returns the computed boundary and interior knots. If these are already both computed (i.e. not nothing), then returns the passed values.