Drake
Drake C++ Documentation
Loading...
Searching...
No Matches
BsplineBasis< T > Class Template Referencefinal

Detailed Description

template<typename T>
class drake::math::BsplineBasis< T >

Given a set of non-descending breakpoints t₀ ≤ t₁ ≤ ⋅⋅⋅ ≤ tₘ, a B-spline basis of order k is a set of n + 1 (where n = m - k) piecewise polynomials of degree k - 1 defined over those breakpoints.

The elements of this set are called "B-splines". The vector (t₀, t₁, ..., tₘ)' is referred to as the "knot vector" of the basis and its elements are referred to as "knots".

At a breakpoint with multiplicity p (i.e. a breakpoint that appears p times in the knot vector), B-splines are guaranteed to have Cᵏ⁻ᵖ⁻¹ continuity.

A B-spline curve using a B-spline basis B, is a parametric curve mapping parameter values in [tₖ₋₁, tₙ₊₁] to a vector space V. For t ∈ [tₖ₋₁, tₙ₊₁] the value of the curve is given by the linear combination of n + 1 control points, pᵢ ∈ V, with the elements of B evaluated at t.

For more information on B-splines and their uses, see (for example) Patrikalakis et al. [1].

[1] https://web.mit.edu/hyperbook/Patrikalakis-Maekawa-Cho/node15.html

#include <drake/math/bspline_basis.h>

Public Member Functions

 BsplineBasis ()
 BsplineBasis (int order, std::vector< T > knots)
 Constructs a B-spline basis with the specified order and knots.
 BsplineBasis (int order, int num_basis_functions, KnotVectorType type=KnotVectorType::kClampedUniform, const T &initial_parameter_value=0, const T &final_parameter_value=1)
 Constructs a B-spline basis with the specified order, num_basis_functions, initial_parameter_value, final_parameter_value, and an auto-generated knot vector of the specified type.
 BsplineBasis (const BsplineBasis< double > &other)
 Conversion constructor.
int order () const
 The order of this B-spline basis (k in the class description).
int degree () const
 The degree of the piecewise polynomials comprising this B-spline basis (k - 1 in the class description).
int num_basis_functions () const
 The number of basis functions in this B-spline basis (n + 1 in the class description).
const std::vector< T > & knots () const
 The knot vector of this B-spline basis (the vector (t₀, t₁, ..., tₘ)' in the class description).
const T & initial_parameter_value () const
 The minimum allowable parameter value for B-spline curves using this basis (tₖ₋₁ in the class description).
const T & final_parameter_value () const
 The maximum allowable parameter value for B-spline curves using this basis (tₙ₊₁ in the class description).
int FindContainingInterval (const T &parameter_value) const
 For a parameter_value = t, the interval that contains it is the pair of knot values [tᵢ, tᵢ₊₁] for the greatest i such that tᵢ ≤ t and tᵢ < final_parameter_value().
std::vector< int > ComputeActiveBasisFunctionIndices (const std::array< T, 2 > &parameter_interval) const
 Returns the indices of the basis functions which may evaluate to non-zero values for some parameter value in parameter_interval; all other basis functions are strictly zero over parameter_interval.
std::vector< int > ComputeActiveBasisFunctionIndices (const T &parameter_value) const
 Returns the indices of the basis functions which may evaluate to non-zero values for parameter_value; all other basis functions are strictly zero at this point.
template<typename T_control_point>
T_control_point EvaluateCurve (const std::vector< T_control_point > &control_points, const T &parameter_value) const
 Evaluates the B-spline curve defined by this and control_points at the given parameter_value.
VectorX< T > EvaluateLinearInControlPoints (const T &parameter_value) const
 Returns the vector, M, such that.
EvaluateBasisFunctionI (int i, const T &parameter_value) const
 Returns the value of the i-th basis function evaluated at parameter_value.
boolean< T > operator== (const BsplineBasis &other) const
boolean< T > operator!= (const BsplineBasis &other) const
template<typename Archive>
void Serialize (Archive *a)
 Passes this object to an Archive.
Implements CopyConstructible, CopyAssignable, MoveConstructible, MoveAssignable
 BsplineBasis (const BsplineBasis &)=default
BsplineBasisoperator= (const BsplineBasis &)=default
 BsplineBasis (BsplineBasis &&)=default
BsplineBasisoperator= (BsplineBasis &&)=default

Constructor & Destructor Documentation

◆ BsplineBasis() [1/6]

template<typename T>
BsplineBasis ( const BsplineBasis< T > & )
default

◆ BsplineBasis() [2/6]

template<typename T>
BsplineBasis ( BsplineBasis< T > && )
default

◆ BsplineBasis() [3/6]

template<typename T>
BsplineBasis ( )

◆ BsplineBasis() [4/6]

template<typename T>
BsplineBasis ( int order,
std::vector< T > knots )

Constructs a B-spline basis with the specified order and knots.

Precondition
knots is sorted in non-descending order.
Exceptions
std::exceptionif knots.size() < 2 * order.

◆ BsplineBasis() [5/6]

template<typename T>
BsplineBasis ( int order,
int num_basis_functions,
KnotVectorType type = KnotVectorType::kClampedUniform,
const T & initial_parameter_value = 0,
const T & final_parameter_value = 1 )

Constructs a B-spline basis with the specified order, num_basis_functions, initial_parameter_value, final_parameter_value, and an auto-generated knot vector of the specified type.

Exceptions
std::exceptionif num_basis_functions < order
Precondition
initial_parameter_value ≤ final_parameter_value

◆ BsplineBasis() [6/6]

template<typename T>
BsplineBasis ( const BsplineBasis< double > & other)
explicit

Conversion constructor.

Constructs an instance of BsplineBasis<T> from a double-valued basis.

Member Function Documentation

◆ ComputeActiveBasisFunctionIndices() [1/2]

template<typename T>
std::vector< int > ComputeActiveBasisFunctionIndices ( const std::array< T, 2 > & parameter_interval) const

Returns the indices of the basis functions which may evaluate to non-zero values for some parameter value in parameter_interval; all other basis functions are strictly zero over parameter_interval.

Precondition
parameter_interval[0] ≤ parameter_interval[1]
parameter_interval[0] ≥ initial_parameter_value()
parameter_interval[1] ≤ final_parameter_value()

◆ ComputeActiveBasisFunctionIndices() [2/2]

template<typename T>
std::vector< int > ComputeActiveBasisFunctionIndices ( const T & parameter_value) const

Returns the indices of the basis functions which may evaluate to non-zero values for parameter_value; all other basis functions are strictly zero at this point.

Precondition
parameter_value ≥ initial_parameter_value()
parameter_value ≤ final_parameter_value()

◆ degree()

template<typename T>
int degree ( ) const

The degree of the piecewise polynomials comprising this B-spline basis (k - 1 in the class description).

◆ EvaluateBasisFunctionI()

template<typename T>
T EvaluateBasisFunctionI ( int i,
const T & parameter_value ) const

Returns the value of the i-th basis function evaluated at parameter_value.

◆ EvaluateCurve()

template<typename T>
template<typename T_control_point>
T_control_point EvaluateCurve ( const std::vector< T_control_point > & control_points,
const T & parameter_value ) const

Evaluates the B-spline curve defined by this and control_points at the given parameter_value.

Parameters
control_pointsControl points of the B-spline curve.
parameter_valueParameter value at which to evaluate the B-spline curve defined by this and control_points.
Precondition
control_points.size() == num_basis_functions()
parameter_value ≥ initial_parameter_value()
parameter_value ≤ final_parameter_value()

◆ EvaluateLinearInControlPoints()

template<typename T>
VectorX< T > EvaluateLinearInControlPoints ( const T & parameter_value) const

Returns the vector, M, such that.

EvaluateCurve(control_points, parameter_value) = control_points * M

where T_control_points==VectorX<T> (so control_points is a matrix). This is useful for writing linear constraints on the control points.

Precondition
parameter_value ≥ initial_parameter_value()
parameter_value ≤ final_parameter_value()

◆ final_parameter_value()

template<typename T>
const T & final_parameter_value ( ) const

The maximum allowable parameter value for B-spline curves using this basis (tₙ₊₁ in the class description).

◆ FindContainingInterval()

template<typename T>
int FindContainingInterval ( const T & parameter_value) const

For a parameter_value = t, the interval that contains it is the pair of knot values [tᵢ, tᵢ₊₁] for the greatest i such that tᵢ ≤ t and tᵢ < final_parameter_value().

This function returns that value of i.

Precondition
parameter_value ≥ initial_parameter_value()
parameter_value ≤ final_parameter_value()

◆ initial_parameter_value()

template<typename T>
const T & initial_parameter_value ( ) const

The minimum allowable parameter value for B-spline curves using this basis (tₖ₋₁ in the class description).

◆ knots()

template<typename T>
const std::vector< T > & knots ( ) const

The knot vector of this B-spline basis (the vector (t₀, t₁, ..., tₘ)' in the class description).

◆ num_basis_functions()

template<typename T>
int num_basis_functions ( ) const

The number of basis functions in this B-spline basis (n + 1 in the class description).

◆ operator!=()

template<typename T>
boolean< T > operator!= ( const BsplineBasis< T > & other) const

◆ operator=() [1/2]

template<typename T>
BsplineBasis & operator= ( BsplineBasis< T > && )
default

◆ operator=() [2/2]

template<typename T>
BsplineBasis & operator= ( const BsplineBasis< T > & )
default

◆ operator==()

template<typename T>
boolean< T > operator== ( const BsplineBasis< T > & other) const

◆ order()

template<typename T>
int order ( ) const

The order of this B-spline basis (k in the class description).

◆ Serialize()

template<typename T>
template<typename Archive>
void Serialize ( Archive * a)

Passes this object to an Archive.

Refer to YAML Serialization for background. This method is only available when T = double.


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