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 functions. The introduction is intended for package developers who would like to use splines2 package at C++ level.
Different with the procedure based functions in the R interface, 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>
// include header file from splines2 package
#include <splines2Armadillo.h>
// [[Rcpp::plugins(cpp11)]]
To use Rcpp::sourceCpp()
, one may need to add [[Rcpp::depends()]]
as follows:
// [[Rcpp::depends(RcppArmadillo)]]
// [[Rcpp::depends(splines2)]]
For ease of demonstration, we assume the following using-directives:
using namespace arma
using namespace splines2
The BernsteinPoly
class is implemented for the generalized Bernstein polynomials.
The main non-default constructor is as follows:
(const vec& x,
BernsteinPolyconst unsigned int degree,
const vec& boundary_knots = vec());
In addition, two explicit constructors are provided for BernsteinPoly*
and SplineBase*
, which set x
, degree
, and boundary_knots
from the objects that the pointers point to.
The main methods are
basis()
for basis matrixderivative()
for derivatives of basis functionsintegral()
for integrals of basis functionsThe specific function signatures are as follows:
(const bool complete_basis = true);
mat basis(const unsigned int derivs = 1,
mat derivativeconst bool complete_basis = true);
(const bool complete_basis = true); mat integral
In addition, we may set and get the specifications through the following setter and getter functions, respectively.
// setter functions
* set_x(const vec&);
BernsteinPoly* set_x(const double);
BernsteinPoly* set_degree(const unsigned int);
BernsteinPoly* set_order(const unsigned int);
BernsteinPoly* set_internal_knots(const vec&); // placeholder, does nothing
BernsteinPoly* set_boundary_knots(const vec&);
BernsteinPoly
// getter functions
();
vec get_xunsigned int get_degree();
unsigned int get_order();
(); vec get_boundary_knots
The setter function returns a pointer to the current object.
A virtual base class named SplineBase
is implemented to support a variety of classes for spline basis functions including BSpline
, MSpline
, ISpline
, CSpline
, NaturalSpline
, PeriodicMSpline
. Their names are self-explanatory.
BSpline
, MSpline
, ISpline
, and CSpline
The BSpline
, MSpline
, ISpline
, and CSpline
classes share the same constructors inherited from the SplineBase
class. There are four constructors in addition to the default constructor.
The first non-default constructor is called when internal knots are explicitly specified. For example,
// 1. specified internal_knots
(const vec& x,
BSplineconst 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 functions (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
. For example,
// 2. specified spline degree of freedom (df)
(const vec& x,
BSplineconst unsigned int spline_df,
const unsigned int degree = 3,
const vec& boundary_knots = vec());
The third non-default constructor is intended for basis functions with an extended knot sequence. The multiplicities of knots in the sequence can be different (but should not be greater than degree + 1
). For example,
// 3. specified degree and (extended) knot sequence
(const vec& x,
BSplineconst unsigned int degree,
const vec& knot_sequence);
The fourth 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. For example,
// 4. create a new object from a base class pointer
(const SplineBase* pSplineBase); BSpline
NaturalSpline
The NaturalSpline
represents the class for natural cubic splines. Thus, its constructors do not allow specification of degree
. The first non-default constructor is called when internal knots are explicitly specified.
// 1. specified internal_knots
(const vec& x,
NaturalSplineconst vec& internal_knots,
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 functions (different with df
in the R interface) is specified. Then the number of internal knots is computed as spline_df - 2
and the placement of internal knots uses quantiles of specified x
.
// 2. specified spline degree of freedom (df)
(const vec& x,
NaturalSplineconst unsigned int spline_df,
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 (x
, internal_knots
, boundary_knots
, etc.) of an existing object.
// 3. create a new object from a base class pointer
(const SplineBase* pSplineBase); NaturalSpline
PeriodicMSpline
The PeriodicMSpline
class is for constructing the periodic M-splines, which provides the same set of non-default constructors with BSpline
except the constructor for directly specifying the knot sequence.
The main methods are
basis()
for spline basis matrixderivative()
for derivatives of spline basisintegral()
for integrals of spline basis (except for the CSpline
class)The specific function signatures are as follows:
(const bool complete_basis = true);
mat basis(const unsigned int derivs = 1,
mat derivativeconst bool complete_basis = true);
(const bool complete_basis = true); mat integral
Similarly, we may set and get the spline specifications through the following setter and getter functions, respectively.
// setter functions
* set_x(const vec&);
SplineBase* set_x(const double);
SplineBase* set_internal_knots(const vec&);
SplineBase* set_boundary_knots(const vec&);
SplineBase* set_knot_sequence(const vec&);
SplineBase* set_degree(const unsigned int);
SplineBase* set_order(const unsigned int);
SplineBase
// getter functions
();
vec get_x();
vec get_internal_knots();
vec get_boundary_knots();
vec get_knot_sequenceunsigned int get_degree();
unsigned int get_order();
unsigned int get_spline_df();
The setter function returns a pointer to the current object so that the specification can be chained for convenience. For example,
{ arma::regspace(0, 0.1, 1) }; // 0, 0.1, ..., 1
vec x { x, 5 }; // df = 5 (and degree = 3, by default)
BSpline obj // change degree to 2 and get basis
{ obj.set_degree(2)->basis() }; mat basis_mat
There is no available integral()
method for CSpline
and no meaningful degree
related methods for NaturalSpline
.