Drake
Drake C++ Documentation
BezierCurve< T > Class Template Referencefinal

Detailed Description

template<typename T>
class drake::trajectories::BezierCurve< T >

A Bézier curve is defined by a set of control points p₀ through pₙ, where n is called the order of the curve (n = 1 for linear, 2 for quadratic, 3 for cubic, etc.).

The first and last control points are always the endpoints of the curve; however, the intermediate control points (if any) generally do not lie on the curve, but the curve is guaranteed to stay within the convex hull of the control points.

See also BsplineTrajectory. A B-spline can be thought of as a composition of overlapping Bézier curves (where each evaluation only depends on a local subset of the control points). In contrast, evaluating a Bézier curve will use all of the control points.

Template Parameters
TThe scalar type, which must be one of the default scalars.

#include <drake/common/trajectories/bezier_curve.h>

Public Member Functions

 BezierCurve ()
 Default initializer. More...
 
 BezierCurve (double start_time, double end_time, const Eigen::Ref< const MatrixX< T >> &control_points)
 Constructs a Bézier curve over the interval t ∈ [start_time, end_time] with control points defined in the columns of control_points. More...
 
 ~BezierCurve () override
 
int order () const
 Returns the order of the curve (1 for linear, 2 for quadratic, etc.). More...
 
BernsteinBasis (int i, const T &time, std::optional< int > order=std::nullopt) const
 Returns the value of the ith basis function of order (1 for linear, 2 for quadratic, etc) evaluated at time. More...
 
const MatrixX< T > & control_points () const
 Returns a const reference to the control points which define the curve. More...
 
Eigen::SparseMatrix< doubleAsLinearInControlPoints (int derivative_order=1) const
 Supports writing optimizations using the control points as decision variables. More...
 
std::unique_ptr< trajectories::Trajectory< T > > Clone () const override
 
MatrixX< T > value (const T &time) const override
 Evaluates the curve at the given time. More...
 
VectorX< symbolic::ExpressionGetExpression (symbolic::Variable time=symbolic::Variable("t")) const
 Extracts the expanded underlying polynomial expression of this curve in terms of variable time. More...
 
void ElevateOrder ()
 Increases the order of the curve by 1. More...
 
Eigen::Index rows () const override
 
Eigen::Index cols () const override
 
start_time () const override
 
end_time () const override
 
Implements CopyConstructible, CopyAssignable, MoveConstructible, MoveAssignable
 BezierCurve (const BezierCurve &)=default
 
BezierCurveoperator= (const BezierCurve &)=default
 
 BezierCurve (BezierCurve &&)=default
 
BezierCurveoperator= (BezierCurve &&)=default
 
- Public Member Functions inherited from Trajectory< T >
virtual ~Trajectory ()
 
MatrixX< T > vector_values (const std::vector< T > &t) const
 If cols()==1, then evaluates the trajectory at each time t, and returns the results as a Matrix with the ith column corresponding to the ith time. More...
 
MatrixX< T > vector_values (const Eigen::Ref< const VectorX< T >> &t) const
 If cols()==1, then evaluates the trajectory at each time t, and returns the results as a Matrix with the ith column corresponding to the ith time. More...
 
bool has_derivative () const
 Returns true iff the Trajectory provides and implementation for EvalDerivative() and MakeDerivative(). More...
 
MatrixX< T > EvalDerivative (const T &t, int derivative_order=1) const
 Evaluates the derivative of this at the given time t. More...
 
std::unique_ptr< Trajectory< T > > MakeDerivative (int derivative_order=1) const
 Takes the derivative of this Trajectory. More...
 

Additional Inherited Members

- Protected Member Functions inherited from Trajectory< T >
 Trajectory ()=default
 
 Trajectory (const Trajectory &)=default
 
Trajectoryoperator= (const Trajectory &)=default
 
 Trajectory (Trajectory &&)=default
 
Trajectoryoperator= (Trajectory &&)=default
 

Constructor & Destructor Documentation

◆ BezierCurve() [1/4]

BezierCurve ( const BezierCurve< T > &  )
default

◆ BezierCurve() [2/4]

BezierCurve ( BezierCurve< T > &&  )
default

◆ BezierCurve() [3/4]

Default initializer.

Constructs an empty Bézier curve over the interval t ∈ [0, 1].

◆ BezierCurve() [4/4]

BezierCurve ( double  start_time,
double  end_time,
const Eigen::Ref< const MatrixX< T >> &  control_points 
)

Constructs a Bézier curve over the interval t ∈ [start_time, end_time] with control points defined in the columns of control_points.

Precondition
end_time >= start_time.

◆ ~BezierCurve()

~BezierCurve ( )
override

Member Function Documentation

◆ AsLinearInControlPoints()

Eigen::SparseMatrix<double> AsLinearInControlPoints ( int  derivative_order = 1) const

Supports writing optimizations using the control points as decision variables.

This method returns the matrix, M, defining the control points of the order derivative in the form:

derivative.control_points() = this.control_points() * M

For instance, since we have

derivative.control_points().col(k) = this.control_points() * M.col(k),

constraining the kth control point of the nth derivative to be in [ub, lb], could be done with:

auto M = curve.AsLinearInControlPoints(n);
for (int i=0; i<curve.rows(); ++i) {
auto c = std::make_shared<solvers::LinearConstraint>(
M.col(k).transpose(), Vector1d(lb(i)), Vector1d(ub(i)));
prog.AddConstraint(c, curve.row(i).transpose());
}

Iterating over the rows of the control points is the natural sparsity pattern here (since M is the same for all rows). For instance, we also have

derivative.control_points().row(k).T = M.T * this.control_points().row(k).T,

or

vec(derivative.control_points().T) = blockMT * vec(this.control_points().T),
blockMT = [ M.T,   0, .... 0 ]
          [   0, M.T, 0, ... ]
          [      ...         ]
          [  ...    , 0, M.T ].
Precondition
derivative_order >= 0.

◆ BernsteinBasis()

T BernsteinBasis ( int  i,
const T &  time,
std::optional< int order = std::nullopt 
) const

Returns the value of the ith basis function of order (1 for linear, 2 for quadratic, etc) evaluated at time.

The default value for the optional argument order is the order() of this.

◆ Clone()

std::unique_ptr<trajectories::Trajectory<T> > Clone ( ) const
overridevirtual
Returns
A deep copy of this Trajectory.

Implements Trajectory< T >.

◆ cols()

Eigen::Index cols ( ) const
overridevirtual
Returns
The number of columns in the matrix returned by value().

Implements Trajectory< T >.

◆ control_points()

const MatrixX<T>& control_points ( ) const

Returns a const reference to the control points which define the curve.

◆ ElevateOrder()

void ElevateOrder ( )

Increases the order of the curve by 1.

A Bézier curve of order n can be converted into a Bézier curve of order n + 1 with the same shape. The control points of this are modified to obtain the equivalent curve.

◆ end_time()

T end_time ( ) const
overridevirtual

Implements Trajectory< T >.

◆ GetExpression()

VectorX<symbolic::Expression> GetExpression ( symbolic::Variable  time = symbolic::Variable("t")) const

Extracts the expanded underlying polynomial expression of this curve in terms of variable time.

◆ operator=() [1/2]

BezierCurve& operator= ( BezierCurve< T > &&  )
default

◆ operator=() [2/2]

BezierCurve& operator= ( const BezierCurve< T > &  )
default

◆ order()

int order ( ) const

Returns the order of the curve (1 for linear, 2 for quadratic, etc.).

◆ rows()

Eigen::Index rows ( ) const
overridevirtual
Returns
The number of rows in the matrix returned by value().

Implements Trajectory< T >.

◆ start_time()

T start_time ( ) const
overridevirtual

Implements Trajectory< T >.

◆ value()

MatrixX<T> value ( const T &  time) const
overridevirtual

Evaluates the curve at the given time.

Warning
If t does not lie in the range [start_time(), end_time()], the trajectory will silently be evaluated at the closest valid value of time to time. For example, value(-1) will return value(0) for a trajectory defined over [0, 1].

Implements Trajectory< T >.


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