ampsci
High-precision calculations for one- and two-valence atomic systems
BSpline

ok

Basis of N B-splines of order K (degree K-1), defined over [0, xmax].

Constructs and evaluates a B-spline basis using the GSL B-spline library (gsl/gsl_bspline.h). See also: Bachau et al., Rep. Prog. Phys. 64, 1815 (2001).

Knot structure:**

  • N - K + 2 breakpoints (internal knots) are placed in [x0, xmax], where x0 is the first non-zero knot. The point 0 is prepended, giving N - K + 2 breakpoints total spanning [0, xmax].
  • The full knot vector repeats the endpoints K times.
  • Breakpoints may be placed on a logarithmic (default), linear, or log-linear scale; see KnotDistro.

    Compatibility:** Supports GSL v2.x (<2.8), and v2.8+. API differences are handled internally via compile-time version macros. Should also support v1.x, but this is difficult to test.

Copy construction and copy assignment are deleted (GSL workspace is not copyable).

#include <BSpline.hpp>

Public Types

enum class  KnotDistro { logarithmic , linear , loglinear }
 Distribution of internal knot points. More...
 

Public Member Functions

 BSpline (std::size_t n, std::size_t k, double x0, double xmax, KnotDistro kd=KnotDistro::logarithmic)
 Construct a basis of N B-splines of order K over [0, xmax].
 
BSplineoperator= (const BSpline &)=delete
 Copy assignment deleted: GSL workspace is not copyable.
 
 BSpline (const BSpline &)=delete
 Copy construction deleted: GSL workspace is not copyable.
 
std::size_t K ()
 Order of the splines, K.
 
std::size_t d ()
 Degree of the splines, d := K-1.
 
std::size_t N ()
 Number of splines, N.
 
std::size_t find_i0 (double x)
 Returns the index i0 of the first non-zero spline at x. Spline indices run [0, N); at most K splines are non-zero at any x, occupying indices [i0, min(i0+K-1, N-1)].
 
std::vector< double > get (double x)
 Returns all N spline values {b_0(x), b_1(x), ..., b_{N-1}(x)}. Returns a vector of zeros if x is outside [0, xmax].
 
std::pair< std::size_t, LinAlg::Matrix< double > > get_nonzero (double x, std::size_t n_deriv)
 Returns the non-zero splines and their derivatives at x.
 
void print_knots () const
 Prints the full knot vector to stdout (N+K values).
 

Member Enumeration Documentation

◆ KnotDistro

enum class BSpline::KnotDistro
strong

Distribution of internal knot points.

Enumerator
logarithmic 

Knots on a logarithmic scale (default; suits bound-state grids)

linear 

Knots evenly spaced.

loglinear 

Log-linear blend; midpoint at xmax/2.

Constructor & Destructor Documentation

◆ BSpline() [1/2]

BSpline::BSpline ( std::size_t  n,
std::size_t  k,
double  x0,
double  xmax,
KnotDistro  kd = KnotDistro::logarithmic 
)
inline

Construct a basis of N B-splines of order K over [0, xmax].

Parameters
nNumber of splines N. Must satisfy N > K.
kOrder K (degree K-1). E.g. K=4 gives cubic splines.
x0First non-zero (internal) knot; must satisfy 0 < x0 < xmax.
xmaxUpper bound of the spline domain.
kdKnot distribution: logarithmic (default), linear, or loglinear.

◆ BSpline() [2/2]

BSpline::BSpline ( const BSpline )
delete

Copy construction deleted: GSL workspace is not copyable.

Member Function Documentation

◆ operator=()

BSpline & BSpline::operator= ( const BSpline )
delete

Copy assignment deleted: GSL workspace is not copyable.

◆ K()

std::size_t BSpline::K ( )
inline

Order of the splines, K.

◆ d()

std::size_t BSpline::d ( )
inline

Degree of the splines, d := K-1.

◆ N()

std::size_t BSpline::N ( )
inline

Number of splines, N.

◆ find_i0()

std::size_t BSpline::find_i0 ( double  x)
inline

Returns the index i0 of the first non-zero spline at x. Spline indices run [0, N); at most K splines are non-zero at any x, occupying indices [i0, min(i0+K-1, N-1)].

◆ get()

std::vector< double > BSpline::get ( double  x)
inline

Returns all N spline values {b_0(x), b_1(x), ..., b_{N-1}(x)}. Returns a vector of zeros if x is outside [0, xmax].

◆ get_nonzero()

std::pair< std::size_t, LinAlg::Matrix< double > > BSpline::get_nonzero ( double  x,
std::size_t  n_deriv 
)
inline

Returns the non-zero splines and their derivatives at x.

Returns {i0, M}, where:

  • i0 is the index of the first non-zero spline,
  • M is a K (n_deriv+1) matrix; M(i, j) is the j-th derivative of spline b[i0 + i] evaluated at x.

At most K splines are non-zero at any x. This is more efficient than get() when only the non-zero splines are needed.

If x > xmax, returns a zero matrix with i0+K > N. If x = xmax, only the last spline is non-zero (lower derivatives may not be zero).

Parameters
xEvaluation point.
n_derivMaximum derivative order to compute (inclusive).
Returns
Pair {i0, M}.

◆ print_knots()

void BSpline::print_knots ( ) const
inline

Prints the full knot vector to stdout (N+K values).


The documentation for this class was generated from the following file: