Loading...
Searching...
No Matches
ndcurves::bezier_curve< Time, Numeric, Safe, Point > Struct Template Reference

#include <ndcurves/bezier_curve.h>

Inheritance diagram for ndcurves::bezier_curve< Time, Numeric, Safe, Point >:
Collaboration diagram for ndcurves::bezier_curve< Time, Numeric, Safe, Point >:

Public Types

typedef Point point_t
 
typedef Eigen::Matrix< Numeric, Eigen::Dynamic, 1 > vector_x_t
 
typedef Eigen::Ref< const vector_x_tvector_x_ref_t
 
typedef Time time_t
 
typedef Numeric num_t
 
typedef curve_constraints< point_tcurve_constraints_t
 
typedef std::vector< point_t, Eigen::aligned_allocator< point_t > > t_point_t
 
typedef t_point_t::const_iterator cit_point_t
 
typedef bezier_curve< Time, Numeric, Safe, Point > bezier_curve_t
 
typedef std::shared_ptr< bezier_curve_tbezier_curve_ptr_t
 
typedef piecewise_curve< Time, Numeric, Safe, point_t, point_t, bezier_curve_tpiecewise_curve_t
 
typedef curve_abc< Time, Numeric, Safe, point_tcurve_abc_t
 
typedef curve_abc_t::curve_ptr_t curve_ptr_t
 
- Public Types inherited from ndcurves::curve_abc< Time, Numeric, Safe, Point, Point_derivate >
typedef Point point_t
 
typedef Point_derivate point_derivate_t
 
typedef Time time_t
 
typedef Numeric num_t
 
typedef curve_abc< Time, Numeric, Safe, point_t, point_derivate_tcurve_t
 
typedef curve_abc< Time, Numeric, Safe, point_derivate_tcurve_derivate_t
 
typedef std::shared_ptr< curve_tcurve_ptr_t
 

Public Member Functions

 bezier_curve ()
 Empty constructor. Curve obtained this way can not perform other class functions.
 
template<typename In >
 bezier_curve (In PointsBegin, In PointsEnd, const time_t T_min=0., const time_t T_max=1., const time_t mult_T=1.)
 Constructor. Given the first and last point of a control points set, create the bezier curve.
 
template<typename In >
 bezier_curve (In PointsBegin, In PointsEnd, const curve_constraints_t &constraints, const time_t T_min=0., const time_t T_max=1., const time_t mult_T=1.)
 Constructor with constraints. This constructor will add 4 points (2 after the first one, 2 before the last one) to ensure that velocity and acceleration constraints are respected.
 
 bezier_curve (const bezier_curve &other)
 
virtual ~bezier_curve ()
 Destructor.
 
virtual point_t operator() (const time_t t) const
 Evaluation of the bezier curve at time t.
 
bool isApprox (const bezier_curve_t &other, const Numeric prec=Eigen::NumTraits< Numeric >::dummy_precision()) const
 isApprox check if other and *this are approximately equals. Only two curves of the same class can be approximately equals, for comparison between different type of curves see isEquivalent
 
virtual bool isApprox (const curve_abc_t *other, const Numeric prec=Eigen::NumTraits< Numeric >::dummy_precision()) const
 
virtual bool operator== (const bezier_curve_t &other) const
 
virtual bool operator!= (const bezier_curve_t &other) const
 
bezier_curve_t compute_derivate (const std::size_t order) const
 Compute the derived curve at order N. Computes the derivative order N, \(\frac{d^Nx(t)}{dt^N}\) of bezier curve of parametric equation x(t).
 
bezier_curve_tcompute_derivate_ptr (const std::size_t order) const
 Compute the derived curve at order N.
 
bezier_curve_t compute_primitive (const std::size_t order, const point_t &init) const
 Compute the primitive of the curve at order N. Computes the primitive at order N of bezier curve of parametric equation \(x(t)\).
At order \(N=1\), the primitve \(X(t)\) of \(x(t)\) is such as \(\frac{dX(t)}{dt} = x(t)\).
 
bezier_curve_t compute_primitive (const std::size_t order) const
 
bezier_curve_tcompute_primitive_ptr (const std::size_t order, const point_t &init) const
 
bezier_curve_t elevate (const std::size_t order) const
 Computes a Bezier curve of order degrees higher than the current curve, but strictly equivalent. Order elevation is required for addition / substraction and other comparison operations.
 
void elevate_self (const std::size_t order)
 Elevate the Bezier curve of order degrees higher than the current curve, but strictly equivalent. Order elevation is required for addition / substraction and other comparison operations.
 
virtual point_t derivate (const time_t t, const std::size_t order) const
 Evaluate the derivative order N of curve at time t. If derivative is to be evaluated several times, it is rather recommended to compute derived curve using compute_derivate.
 
point_t evalBernstein (const Numeric t) const
 Evaluate all Bernstein polynomes for a certain degree. A bezier curve with N control points is represented by : \(x(t) = \sum_{i=0}^{N} B_i^N(t) P_i\) with \( B_i^N(t) = \binom{N}{i}t^i (1-t)^{N-i} \).
Warning: the horner scheme is about 100 times faster than this method.
This method will probably be removed in the future as the computation of bernstein polynomial is very costly.
 
- Public Member Functions inherited from ndcurves::curve_abc< Time, Numeric, Safe, Point, Point_derivate >
 curve_abc ()
 Constructor.
 
virtual ~curve_abc ()
 Destructor.
 
bool isEquivalent (const curve_t *other, const Numeric prec=Eigen::NumTraits< Numeric >::dummy_precision(), const size_t order=5) const
 isEquivalent check if other and *this are approximately equal by values, given a precision threshold. This test is done by discretizing both curves and evaluating them and their derivatives.
 
virtual bool isApprox (const curve_t *other, const Numeric prec=Eigen::NumTraits< Numeric >::dummy_precision()) const =0
 isApprox check if other and *this are approximately equal given a precision threshold Only two curves of the same class can be approximately equal, for comparison between different type of curves see isEquivalent.
 
virtual std::size_t dim () const =0
 Get dimension of curve.
 
virtual time_t min () const =0
 Get the minimum time for which the curve is defined.
 
virtual time_t max () const =0
 Get the maximum time for which the curve is defined.
 
virtual std::size_t degree () const =0
 Get the degree of the curve.
 
std::pair< time_t, time_ttimeRange ()
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
- Public Member Functions inherited from ndcurves::serialization::Serializable
template<class Derived >
void loadFromText (const std::string &filename)
 Loads a Derived object from a text file.
 
template<class Derived >
void saveAsText (const std::string &filename) const
 Saved a Derived object as a text file.
 
template<class Derived >
void loadFromXML (const std::string &filename, const std::string &tag_name)
 Loads a Derived object from an XML file.
 
template<class Derived >
void saveAsXML (const std::string &filename, const std::string &tag_name) const
 Saved a Derived object as an XML file.
 
template<class Derived >
void loadFromBinary (const std::string &filename)
 Loads a Derived object from an binary file.
 
template<class Derived >
void saveAsBinary (const std::string &filename) const
 Saved a Derived object as an binary file.
 

Member Typedef Documentation

◆ bezier_curve_ptr_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef std::shared_ptr<bezier_curve_t> ndcurves::bezier_curve< Time, Numeric, Safe, Point >::bezier_curve_ptr_t

◆ bezier_curve_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef bezier_curve<Time, Numeric, Safe, Point> ndcurves::bezier_curve< Time, Numeric, Safe, Point >::bezier_curve_t

◆ cit_point_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef t_point_t::const_iterator ndcurves::bezier_curve< Time, Numeric, Safe, Point >::cit_point_t

◆ curve_abc_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef curve_abc<Time, Numeric, Safe, point_t> ndcurves::bezier_curve< Time, Numeric, Safe, Point >::curve_abc_t

◆ curve_constraints_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef curve_constraints<point_t> ndcurves::bezier_curve< Time, Numeric, Safe, Point >::curve_constraints_t

◆ curve_ptr_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef curve_abc_t::curve_ptr_t ndcurves::bezier_curve< Time, Numeric, Safe, Point >::curve_ptr_t

◆ num_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef Numeric ndcurves::bezier_curve< Time, Numeric, Safe, Point >::num_t

◆ piecewise_curve_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef piecewise_curve<Time, Numeric, Safe, point_t, point_t, bezier_curve_t> ndcurves::bezier_curve< Time, Numeric, Safe, Point >::piecewise_curve_t

◆ point_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef Point ndcurves::bezier_curve< Time, Numeric, Safe, Point >::point_t

◆ t_point_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef std::vector<point_t, Eigen::aligned_allocator<point_t> > ndcurves::bezier_curve< Time, Numeric, Safe, Point >::t_point_t

◆ time_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef Time ndcurves::bezier_curve< Time, Numeric, Safe, Point >::time_t

◆ vector_x_ref_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef Eigen::Ref<const vector_x_t> ndcurves::bezier_curve< Time, Numeric, Safe, Point >::vector_x_ref_t

◆ vector_x_t

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
typedef Eigen::Matrix<Numeric, Eigen::Dynamic, 1> ndcurves::bezier_curve< Time, Numeric, Safe, Point >::vector_x_t

Constructor & Destructor Documentation

◆ bezier_curve() [1/4]

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
ndcurves::bezier_curve< Time, Numeric, Safe, Point >::bezier_curve ( )
inline

Empty constructor. Curve obtained this way can not perform other class functions.

◆ bezier_curve() [2/4]

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
template<typename In >
ndcurves::bezier_curve< Time, Numeric, Safe, Point >::bezier_curve ( In  PointsBegin,
In  PointsEnd,
const time_t  T_min = 0.,
const time_t  T_max = 1.,
const time_t  mult_T = 1. 
)
inline

Constructor. Given the first and last point of a control points set, create the bezier curve.

Parameters
PointsBegin: an iterator pointing to the first element of a control point container.
PointsEnd: an iterator pointing to the last element of a control point container.
T_min: lower bound of time, curve will be defined for time in [T_min, T_max].
T_max: upper bound of time, curve will be defined for time in [T_min, T_max].
mult_T: ... (default value is 1.0).

◆ bezier_curve() [3/4]

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
template<typename In >
ndcurves::bezier_curve< Time, Numeric, Safe, Point >::bezier_curve ( In  PointsBegin,
In  PointsEnd,
const curve_constraints_t constraints,
const time_t  T_min = 0.,
const time_t  T_max = 1.,
const time_t  mult_T = 1. 
)
inline

Constructor with constraints. This constructor will add 4 points (2 after the first one, 2 before the last one) to ensure that velocity and acceleration constraints are respected.

Parameters
PointsBegin: an iterator pointing to the first element of a control point container.
PointsEnd: an iterator pointing to the last element of a control point container.
constraints: constraints applying on start / end velocities and acceleration.
T_min: lower bound of time, curve will be defined for time in [T_min, T_max].
T_max: upper bound of time, curve will be defined for time in [T_min, T_max].
mult_T: ... (default value is 1.0).

◆ bezier_curve() [4/4]

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
ndcurves::bezier_curve< Time, Numeric, Safe, Point >::bezier_curve ( const bezier_curve< Time, Numeric, Safe, Point > &  other)
inline

◆ ~bezier_curve()

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
virtual ndcurves::bezier_curve< Time, Numeric, Safe, Point >::~bezier_curve ( )
inlinevirtual

Destructor.

Member Function Documentation

◆ compute_derivate()

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
bezier_curve_t ndcurves::bezier_curve< Time, Numeric, Safe, Point >::compute_derivate ( const std::size_t  order) const
inline

Compute the derived curve at order N. Computes the derivative order N, \(\frac{d^Nx(t)}{dt^N}\) of bezier curve of parametric equation x(t).

Parameters
order: order of derivative.
Returns
\(\frac{d^Nx(t)}{dt^N}\) derivative order N of the curve.

◆ compute_derivate_ptr()

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
bezier_curve_t * ndcurves::bezier_curve< Time, Numeric, Safe, Point >::compute_derivate_ptr ( const std::size_t  order) const
inlinevirtual

Compute the derived curve at order N.

Parameters
order: order of derivative.
Returns
A pointer to \(\frac{d^Nx(t)}{dt^N}\) derivative order N of the curve.

Implements ndcurves::curve_abc< Time, Numeric, Safe, Point, Point_derivate >.

◆ compute_primitive() [1/2]

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
bezier_curve_t ndcurves::bezier_curve< Time, Numeric, Safe, Point >::compute_primitive ( const std::size_t  order) const
inline

◆ compute_primitive() [2/2]

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
bezier_curve_t ndcurves::bezier_curve< Time, Numeric, Safe, Point >::compute_primitive ( const std::size_t  order,
const point_t init 
) const
inline

Compute the primitive of the curve at order N. Computes the primitive at order N of bezier curve of parametric equation \(x(t)\).
At order \(N=1\), the primitve \(X(t)\) of \(x(t)\) is such as \(\frac{dX(t)}{dt} = x(t)\).

Parameters
order: order of the primitive.
init: constant valuefor the first point of the primitive (can tipycally be zero)
Returns
primitive at order N of x(t).

◆ compute_primitive_ptr()

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
bezier_curve_t * ndcurves::bezier_curve< Time, Numeric, Safe, Point >::compute_primitive_ptr ( const std::size_t  order,
const point_t init 
) const
inline

◆ derivate()

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
virtual point_t ndcurves::bezier_curve< Time, Numeric, Safe, Point >::derivate ( const time_t  t,
const std::size_t  order 
) const
inlinevirtual

Evaluate the derivative order N of curve at time t. If derivative is to be evaluated several times, it is rather recommended to compute derived curve using compute_derivate.

Parameters
order: order of derivative.
t: time when to evaluate the curve.
Returns
\(\frac{d^Nx(t)}{dt^N}\) point corresponding on derived curve of order N at time t.

Implements ndcurves::curve_abc< Time, Numeric, Safe, Point, Point_derivate >.

◆ elevate()

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
bezier_curve_t ndcurves::bezier_curve< Time, Numeric, Safe, Point >::elevate ( const std::size_t  order) const
inline

Computes a Bezier curve of order degrees higher than the current curve, but strictly equivalent. Order elevation is required for addition / substraction and other comparison operations.

Parameters
order: number of order the curve must be updated
Returns
An equivalent Bezier, with one more degree.

◆ elevate_self()

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
void ndcurves::bezier_curve< Time, Numeric, Safe, Point >::elevate_self ( const std::size_t  order)
inline

Elevate the Bezier curve of order degrees higher than the current curve, but strictly equivalent. Order elevation is required for addition / substraction and other comparison operations.

Parameters
order: number of order the curve must be updated

◆ evalBernstein()

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
point_t ndcurves::bezier_curve< Time, Numeric, Safe, Point >::evalBernstein ( const Numeric  t) const
inline

Evaluate all Bernstein polynomes for a certain degree. A bezier curve with N control points is represented by : \(x(t) = \sum_{i=0}^{N} B_i^N(t) P_i\) with \( B_i^N(t) = \binom{N}{i}t^i (1-t)^{N-i} \).
Warning: the horner scheme is about 100 times faster than this method.
This method will probably be removed in the future as the computation of bernstein polynomial is very costly.

Parameters
t: time when to evaluate the curve.
Returns
\(x(t)\) point corresponding on curve at time t.

◆ isApprox() [1/2]

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
bool ndcurves::bezier_curve< Time, Numeric, Safe, Point >::isApprox ( const bezier_curve_t other,
const Numeric  prec = Eigen::NumTraits<Numeric>::dummy_precision() 
) const
inline

isApprox check if other and *this are approximately equals. Only two curves of the same class can be approximately equals, for comparison between different type of curves see isEquivalent

Parameters
otherthe other curve to check
precthe precision threshold, default Eigen::NumTraits<Numeric>::dummy_precision()
Returns
true if the two curves are approximately equals

◆ isApprox() [2/2]

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
virtual bool ndcurves::bezier_curve< Time, Numeric, Safe, Point >::isApprox ( const curve_abc_t other,
const Numeric  prec = Eigen::NumTraits<Numeric>::dummy_precision() 
) const
inlinevirtual

◆ operator!=()

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
virtual bool ndcurves::bezier_curve< Time, Numeric, Safe, Point >::operator!= ( const bezier_curve_t other) const
inlinevirtual

◆ operator()()

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
virtual point_t ndcurves::bezier_curve< Time, Numeric, Safe, Point >::operator() ( const time_t  t) const
inlinevirtual

Evaluation of the bezier curve at time t.

Parameters
t: time when to evaluate the curve.
Returns
\(x(t)\) point corresponding on curve at time t.

Implements ndcurves::curve_abc< Time, Numeric, Safe, Point, Point_derivate >.

◆ operator==()

template<typename Time = double, typename Numeric = Time, bool Safe = false, typename Point = Eigen::Matrix<Numeric, Eigen::Dynamic, 1>>
virtual bool ndcurves::bezier_curve< Time, Numeric, Safe, Point >::operator== ( const bezier_curve_t other) const
inlinevirtual

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