In this package vignette, we introduce how to use the C++ header-only library that splines2 contains with the Rcpp package (Eddelbuettel 2013) for constructing spline basis. The introduction is intended for package developers who would like to use splines2 package at C++ level.
Different with the procedure based functions at R level, splines2 provides several spline classes in its C++ interface for ease of usage and maintenance. The implementations use the Armadillo (Sanderson 2016) library with help of RcppArmadillo (Eddelbuettel and Sanderson 2014) and require C++11. We may include the header file named splines2Armadillo.h
to get the access to all the classes and implementations in the name space splines2
.
#include <RcppArmadillo.h>
// [[Rcpp::plugins(cpp11)]]
// include header file from splines2 package
#include <splines2Armadillo.h>
// for ease of demonstration
using arma
using splines2
The BSpline
class is for creating B-spline basis.
There are mainly three constructors in addition to the default constructor: BSpline()
.
The first non-default constructor is called when internal knots are explicitly specified.
// 1. specified internal_knots
const vec& x,
BSpline(const vec& internal_knots,
const unsigned int degree = 3,
const vec& boundary_knots = vec())
The second non-default constructor is called when an unsigned integer representing the degree of freedom of the complete spline basis (different with df
in the R interface) is specified. Then the number of internal knots is computed as spline_df - degree - 1
and the placement of internal knots uses quantiles of specified x
.
// 2. specified spline degree of freedom (df)
const vec& x,
BSpline(const unsigned int spline_df,
const unsigned int degree = 3,
const vec& boundary_knots = vec())
The third non-default constructor is explicit and takes a pointer to a base class object. It can be useful when we want to create a new object using the same specification (degree, internal_knots, boundary_knots, etc.) of an existing object.
// 3. create a new object from a base class pointer
const SplineBase* pSplineBase) BSpline(
The main methods are
basis()
for spline basis matrixderivative()
for derivatives of spline basesintegral()
for integrals of spline basesThe specific function signatures are as follows:
const bool complete_basis = true)
mat basis(const unsigned int derivs = 1,
mat derivative(const bool complete_basis = true)
const bool complete_basis = true) mat integral(
In addition, we may set and get the spline specifications through the following setter and getter functions, respectively.
// setter functions
const vec&)
set_x(const double)
set_x(const vec&)
set_internal_knots(const vec&)
set_boundary_knots(const unsigned int)
set_degree(const unsigned int)
set_order(
// getter functions
get_x()
get_internal_knots()
get_boundary_knots()
get_degree()
get_order() get_spline_df()
The setter function returns a pointer to the current object so that the specification can be chained for convenience. For example,
0, 0.1, 1) }; // 0, 0.1, ..., 1
vec x { arma::regspace(5 }; // degree = 3 dy default
BSpline obj { x, // change degree to 2 and get basis
2)->basis() }; mat basis_mat { obj.set_degree(
The class MSpline
for M-splines, ISpline
for I-splines, and CSpline
for C-splines have the exactly same constructors and function members with BSpline
except there is no publicly available integral()
method for CSpline
.
The BernsteinPoly
class is implemented for the generalized Bernstein polynomials.
The main non-default constructor is as follows:
const vec& x,
BernsteinPoly(const unsigned int degree,
const vec& boundary_knots = vec())
Same with BSpline
, the main methods are
basis()
for basis matrixderivative()
for derivatives of basesintegral()
for integrals of basesThe specific function signatures are as follows:
const bool complete_basis = true)
mat basis(const unsigned int derivs = 1,
mat derivative(const bool complete_basis = true)
const bool complete_basis = true) mat integral(
In addition, we may similarly set and get the specifications through the following setter and getter functions, respectively.
// setter functions
const vec&)
set_x(const double)
set_x(const unsigned int)
set_degree(const unsigned int)
set_order(const vec&)
set_boundary_knots(
// getter functions
get_x()
get_degree()
get_order() get_boundary_knots()
The setter function also returns a pointer to the current object.
Eddelbuettel, Dirk. 2013. Seamless R and C++ Integration with Rcpp. Springer.
Eddelbuettel, Dirk, and Conrad Sanderson. 2014. “RcppArmadillo: Accelerating R with High-Performance C++ Linear Algebra.” Computational Statistics and Data Analysis 71: 1054–63. http://dx.doi.org/10.1016/j.csda.2013.02.005.
Sanderson, Conrad. 2016. “Armadillo: An Open Source C++ Linear Algebra Library for Fast Prototyping and Computationally Intensive Experiments.” Journal of Open Source Software 1: 26.