smooth
A C++ library for Lie theory
Loading...
Searching...
No Matches
Public Member Functions | List of all members
BSpline< K, G > Class Template Reference

Cardinal Bspline on a Lie group. More...

#include <bspline.hpp>

Public Member Functions

 BSpline ()
 Construct a constant bspline defined on [0, 1) equal to identity.
 
 BSpline (double t0, double dt, std::vector< G > &&ctrl_pts)
 Create a BSpline.
 
template<std::ranges::range Rv>
requires (std::is_same_v<std::ranges::range_value_t<Rv>, G>)
 BSpline (double t0, double dt, const Rv &ctrl_pts)
 Create a BSpline.
 
 BSpline (const BSpline &)=default
 Copy constructor.
 
 BSpline (BSpline &&)=default
 Move constructor.
 
BSplineoperator= (const BSpline &)=default
 Copy assignment.
 
BSplineoperator= (BSpline &&)=default
 Move assignment.
 
 ~BSpline ()=default
 Descructor.
 
double dt () const
 Distance between knots.
 
double t_min () const
 Minimal time for which spline is defined.
 
double t_max () const
 Maximal time for which spline is defined.
 
const std::vector< G > & ctrl_pts () const
 Access spline control points.
 
template<typename S = double>
CastT< S, G > operator() (const S &t, OptTangent< CastT< S, G > > vel={}, OptTangent< CastT< S, G > > acc={}) const
 Evaluate BSpline.
 

Detailed Description

template<int K, LieGroup G>
class BSpline< K, G >

Cardinal Bspline on a Lie group.

The curve is defined by

\[ g(t) = g_0 * \exp(\tilde B_1(t) v_1) * ... \exp(\tilde B_N(t) v_N) \]

where \(\tilde B_i(t)\) are cumulative Bspline basis functions and \(v_i = g_i \ominus g_{i-1}\) are the control point differences. The control points - knot time correspondence is as follows

KNOT  -K  -K+1   -K+2  ...    0    1   ...  N-K
CTRL   0     1      2  ...    K  K+1          N
                              ^               ^
                            t_min           t_max

The first K ctrl_pts are exterior points and are outside the support of the spline, which means that the spline is defined on \( [t_{min}, t_{max}] = [t0, (N-K)*dt] \).

For interpolation purposes use an odd spline degree and set

t0 = (timestamp of first control point) + dt*(K-1)/2

which aligns control points with the maximum of the corresponding basis function.

Definition at line 47 of file bspline.hpp.

Constructor & Destructor Documentation

◆ BSpline() [1/2]

template<int K, LieGroup G>
BSpline< K, G >::BSpline ( double  t0,
double  dt,
std::vector< G > &&  ctrl_pts 
)

Create a BSpline.

Parameters
t0start of spline
dtdistance between spline knots
ctrl_ptsspline control points

◆ BSpline() [2/2]

template<int K, LieGroup G>
template<std::ranges::range Rv>
requires (std::is_same_v<std::ranges::range_value_t<Rv>, G>)
BSpline< K, G >::BSpline ( double  t0,
double  dt,
const Rv &  ctrl_pts 
)

Create a BSpline.

Template Parameters
Rrange type
Parameters
t0start of spline
dtdistance between spline knots
ctrl_ptsspline control points

Member Function Documentation

◆ operator()()

template<int K, LieGroup G>
template<typename S = double>
CastT< S, G > BSpline< K, G >::operator() ( const S &  t,
OptTangent< CastT< S, G > >  vel = {},
OptTangent< CastT< S, G > >  acc = {} 
) const

Evaluate BSpline.

Template Parameters
Stime type
Parameters
[in]ttime point to evaluate at
[out]veloutput body velocity at evaluation time
[out]accoutput body acceleration at evaluation time
Returns
spline value at time t
Note
Input t is clamped to spline interval of definition

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