Using splines2 with Rcpp

Wenjie Wang

2022-08-12

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 by adding splines2 to the LinkingTo field of package DESCRIPTION file.

Header File and Name Space

Different with the procedure-based functions in the R interface, the C++ interface follows the commonly-used object-oriented design in C++ 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 assume that is enabled for compilation henceforth. 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

Generalized Bernstein Polynomials

The BernsteinPoly class is implemented for the generalized Bernstein polynomials.

Constructors

The main non-default constructor is as follows:

BernsteinPoly(const vec& x,
              const 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.

Function Members

The main methods are

The specific function signatures are as follows:

mat basis(const bool complete_basis = true);
mat derivative(const unsigned int derivs = 1,
               const bool complete_basis = true);
mat integral(const bool complete_basis = true);

In addition, we may set and get the specifications through the following setter and getter functions, respectively.

// setter functions
BernsteinPoly* 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&);

// getter functions
vec get_x();
unsigned int get_degree();
unsigned int get_order();
vec get_boundary_knots();

The setter function returns a pointer to the current object.

Classes for Spline Basis Functions

A virtual base class named SplineBase is implemented to support a variety of classes for spline basis functions including

Constructors of 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 invoked when internal knots are explicitly specified as the second argument. Taking B-splines as an example, the first non-default constructor of a BSpline object is

// 1. specified internal_knots
BSpline(const vec& x,
        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 is specified as the second argument, which represents 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 within boundary.

// 2. specified spline degree of freedom (df)
BSpline(const vec& x,
        const 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, where the multiplicities of the knots can be more than one.

// 3. specified degree and (extended) knot sequence
BSpline(const vec& x,
        const unsigned int degree,
        const vec& knot_sequence);

The fourth non-default constructor is explicit and takes a pointer to a base class object, which can be useful when we want to create a new object using the same specification (x, degree, internal_knots, and boundary_knots) of an existing object (not necessarily a BSpline object).

// 4. create a new object from a base class pointer
BSpline(const SplineBase* pSplineBase);

This constructor also allows us to easily switch between different types of splines. For example, we create a BSpline object bsp_obj form an existing MSpline object msp_obj with the same specification as follows:

BSpline bsp_obj { &msp_obj };

Constructors of 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
NaturalSpline(const vec& x,
              const 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)
NaturalSpline(const vec& x,
              const 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
NaturalSpline(const SplineBase* pSplineBase);

Constructors of 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.

Function Members

The main methods are

The specific function signatures are as follows:

mat basis(const bool complete_basis = true);
mat derivative(const unsigned int derivs = 1,
               const bool complete_basis = true);
mat integral(const bool complete_basis = true);

Similarly, we may set and get the spline specifications through the following setter and getter functions, respectively.

// setter functions
SplineBase* 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);

// getter functions
vec get_x();
vec get_internal_knots();
vec get_boundary_knots();
vec get_knot_sequence();
unsigned 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,

vec x { arma::regspace(0, 0.1, 1) }; // 0, 0.1, ..., 1
BSpline obj { x, 5 };                // df = 5 (and degree = 3, by default)
// change degree to 2 and get basis
mat basis_mat { obj.set_degree(2)->basis() };

The corresponding first derivatives and integrals of the basis functions can be obtained as follows:

mat derivative_mat { bs.derivative() };
mat integral_mat { bs.integral() };

Notice that there is no available integral() method for CSpline and no meaningful degree related methods for NaturalSpline.

Reference

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.