BaseHermite Class ReferenceΒΆ

Splines: BaseHermite Class Reference
Splines
BaseHermite Class Reference
Inheritance diagram for BaseHermite:

Public Member Functions

function BaseHermite ()
 
function base (in ignoredArg, in varargin)
 
function base_D (in ignoredArg, in varargin)
 
function base_DD (in ignoredArg, in varargin)
 
function base_DDD (in ignoredArg, in varargin)
 
function eval (in ignoredArg, in varargin)
 
function eval_D (in ignoredArg, in varargin)
 
function eval_DD (in ignoredArg, in varargin)
 
function eval_DDD (in ignoredArg, in varargin)
 
function base5 (in ignoredArg, in varargin)
 
function base5_D (in ignoredArg, in varargin)
 
function base5_DD (in ignoredArg, in varargin)
 
function base5_DDD (in ignoredArg, in varargin)
 
function base5_DDDD (in ignoredArg, in varargin)
 
function base5_DDDDD (in ignoredArg, in varargin)
 
function eval5 (in ignoredArg, in varargin)
 
function eval5_D (in ignoredArg, in varargin)
 
function eval5_DD (in ignoredArg, in varargin)
 
function eval5_DDD (in ignoredArg, in varargin)
 
function eval5_DDDD (in ignoredArg, in varargin)
 
function eval5_DDDDD (in ignoredArg, in varargin)
 
function hermite_to_bezier (in ignoredArg, in p0, in p1, in t0, in t1)
 
function bezier_to_hermite (in ignoredArg, in p0, in p1, in p2, in p3)
 
function L2_first_derivative (in ignoredArg)
 
function L2_second_derivative (in ignoredArg)
 
function L2_third_derivative (in ignoredArg)
 
function approximate_length (in ignoredArg, in varargin)
 
function cut (in ignoredArg, in varargin)
 

Constructor & Destructor Documentation

◆ BaseHermite()

function BaseHermite::BaseHermite ( )

Build a matlab object storing Hermite base evaluator

Member Function Documentation

◆ approximate_length()

function BaseHermite::approximate_length ( in ignoredArg,
in varargin )

Approximate the length of the cubic polynomial \( \mathbf{p}(t) \) defined on Hermite data:

values = object.approximate_length( t, P0, P1, T0, T1 )
values = object.approximate_length( t, P0, P1, T0, T1, H )

The length is approximated usin 100 linear segment.

◆ base()

function BaseHermite::base ( in ignoredArg,
in varargin )

Evaluate an Hermite base (cubic degree) at point(s) t for the interval \( [0,1] \) (or \( [0,H] \) if second argument is present).

BASE = object.base( t ) % base for the interval [0,1]
BASE = object.base( t, H ) % base for the interval [0,H]

BASE is a matrix length(t) x 4 whose columns are the values of the Hermite base.

\begin{eqnarray*} h_1(t) &=& x^2(3-2x) \\ h_0(t) &=& 1-h_1(t) \\ h_2(t) &=& x(x(x-2)+1) \\ h_3(t) &=& x^2(x-1) \end{eqnarray*}

basis can be returned in separated vectors

[h0,h1,h2,h3] = object.base( t ) % base for the interval [0,1]
[h0,h1,h2,h3] = object.base( t, H ) % base for the interval [0,H]

◆ base5()

function BaseHermite::base5 ( in ignoredArg,
in varargin )

Evaluate an Hermite base (quintic degree) at point(s) t for the interval \( [0,1] \) (or \( [0,H] \) if second argument is present).

BASE = object.base5( t ) % base for the interval [0,1]
BASE = object.base5( t, H ) % base for the interval [0,H]

BASE is a matrix length(t) x 6 whose columns are the values of the Hermite base that are 6 polynomials of degree 5 with the properties

\[ \begin{array}{cccccc} h_1(0) = 1 & h_1(1) = 0 & h'_1(0) = 1 & h'_1(1) = 0 & h''_1(1) = 1 & h''_1(1) = 0 \\ h_0(0) = 0 & h_0(1) = 1 & h'_0(0) = 0 & h'_0(1) = 0 & h''_0(1) = 0 & h''_0(1) = 0 \\ h_2(0) = 0 & h_2(1) = 0 & h'_2(0) = 1 & h'_2(1) = 0 & h''_2(1) = 0 & h''_2(1) = 0 \\ h_3(0) = 0 & h_3(1) = 0 & h'_3(0) = 0 & h'_3(1) = 1 & h''_3(1) = 0 & h''_3(1) = 0 \\ h_4(0) = 0 & h_4(1) = 0 & h'_4(0) = 0 & h'_4(1) = 0 & h''_4(1) = 1 & h''_4(1) = 0 \\ h_5(0) = 0 & h_5(1) = 0 & h'_5(0) = 0 & h'_5(1) = 0 & h''_5(1) = 0 & h''_5(1) = 1 \\ \end{array} \]

basis can be returned in separated vectors

[h0,h1,h2,h3,h4,h5] = object.base5( t ) % base for the interval [0,1]
[h0,h1,h2,h3,h4,h5] = object.base5( t, H ) % base for the interval [0,H]

◆ base5_D()

function BaseHermite::base5_D ( in ignoredArg,
in varargin )

Evaluate an Hermite base derivatives (quintic degree) at point(s) t for the interval \( [0,1] \) (or \( [0,H] \) if second argument is present).

BASE = object.base5_D( t ) % base for the interval [0,1]
BASE = object.base5_D( t, H ) % base for the interval [0,H]

BASE is a matrix length(t) x 6 whose columns are the values of the Hermite base.

basis can be returned in separated vectors

[h0,h1,h2,h3,h4,h5] = object.base5_D( t ) % base for the interval [0,1]
[h0,h1,h2,h3,h4,h5] = object.base5_D( t, H ) % base for the interval [0,H]

◆ base5_DD()

function BaseHermite::base5_DD ( in ignoredArg,
in varargin )

Evaluate an Hermite base second derivatives (quintic degree) at point(s) t for the interval \( [0,1] \) (or \( [0,H] \) if second argument is present).

BASE = object.base5_DD( t ) % base for the interval [0,1]
BASE = object.base5_DD( t, H ) % base for the interval [0,H]

BASE is a matrix length(t) x 6 whose columns are the values of the Hermite base.

basis can be returned in separated vectors

[h0,h1,h2,h3,h4,h5] = object.base5_DD( t ) % base for the interval [0,1]
[h0,h1,h2,h3,h4,h5] = object.base5_DD( t, H ) % base for the interval [0,H]

◆ base5_DDD()

function BaseHermite::base5_DDD ( in ignoredArg,
in varargin )

Evaluate an Hermite base third derivatives (quintic degree) at point(s) t for the interval \( [0,1] \) (or \( [0,H] \) if second argument is present).

BASE = object.base5_DDD( t ) % base for the interval [0,1]
BASE = object.base5_DDD( t, H ) % base for the interval [0,H]

BASE is a matrix length(t) x 6 whose columns are the values of the Hermite base.

basis can be returned in separated vectors

[h0,h1,h2,h3,h4,h5] = object.base5_DDD( t ) % base for the interval [0,1]
[h0,h1,h2,h3,h4,h5] = object.base5_DDD( t, H ) % base for the interval [0,H]

◆ base5_DDDD()

function BaseHermite::base5_DDDD ( in ignoredArg,
in varargin )

Evaluate an Hermite base 4th derivatives (quintic degree) at point(s) t for the interval \( [0,1] \) (or \( [0,H] \) if second argument is present).

BASE = object.base5_DDDD( t ) % base for the interval [0,1]
BASE = object.base5_DDDD( t, H ) % base for the interval [0,H]

BASE is a matrix length(t) x 6 whose columns are the values of the Hermite base.

basis can be returned in separated vectors

[h0,h1,h2,h3,h4,h5] = object.base5_DDDD( t ) % base for the interval [0,1]
[h0,h1,h2,h3,h4,h5] = object.base5_DDDD( t, H ) % base for the interval [0,H]

◆ base5_DDDDD()

function BaseHermite::base5_DDDDD ( in ignoredArg,
in varargin )

Evaluate an Hermite base 5th derivatives (quintic degree) at point(s) t for the interval \( [0,1] \) (or \( [0,H] \) if second argument is present).

BASE = object.base5_DDDDD( t ) % base for the interval [0,1]
BASE = object.base5_DDDDD( t, H ) % base for the interval [0,H]

BASE is a matrix length(t) x 6 whose columns are the values of the Hermite base.

basis can be returned in separated vectors

[h0,h1,h2,h3,h4,h5] = object.base5_DDDDD( t ) % base for the interval [0,1]
[h0,h1,h2,h3,h4,h5] = object.base5_DDDDD( t, H ) % base for the interval [0,H]

◆ base_D()

function BaseHermite::base_D ( in ignoredArg,
in varargin )

Evaluate an Hermite base derivative at point(s) t for the interval \( [0,1] \) (or \( [0,H] \) if second argument is present).

BASE = object.base_D( t ) % base derivative for the interval [0,1]
BASE = object.base_D( t, H ) % base derivative for the interval [0,H]

BASE is a matrix length(t) x 4 whose columns are the values of the Hermite base.

basis can be returned in separated vectors

[h0,h1,h2,h3] = object.base_D( t ) % base derivative for the interval [0,1]
[h0,h1,h2,h3] = object.base_D( t, H ) % base derivative for the interval [0,H]

◆ base_DD()

function BaseHermite::base_DD ( in ignoredArg,
in varargin )

Evaluate an Hermite base second derivative at point(s) t for the interval \( [0,1] \) (or \( [0,H] \) if second argument is present).

BASE = object.base_DD( t ) % base second derivative for the interval [0,1]
BASE = object.base_DD( t, H ) % base second derivative for the interval [0,H]

BASE is a matrix length(t) x 4 whose columns are the values of the Hermite base.

basis can be returned in separated vectors

[h0,h1,h2,h3] = object.base_DD( t ) % base second derivative for the interval [0,1]
[h0,h1,h2,h3] = object.base_DD( t, H ) % base second derivative for the interval [0,H]

◆ base_DDD()

function BaseHermite::base_DDD ( in ignoredArg,
in varargin )

Evaluate an Hermite base third derivative at point(s) t for the interval \( [0,1] \) (or \( [0,H] \) if second argument is present).

BASE = object.base_DDD( t ) % base third derivative for the interval [0,1]
BASE = object.base_DDD( t, H ) % base third derivative for the interval [0,H]

BASE is a matrix length(t) x 4 whose columns are the values of the Hermite base.

basis can be returned in separated vectors

[h0,h1,h2,h3] = object.base_DDD( t ) % base third derivative for the interval [0,1]
[h0,h1,h2,h3] = object.base_DDD( t, H ) % base third derivative for the interval [0,H]

◆ bezier_to_hermite()

function BaseHermite::bezier_to_hermite ( in ignoredArg,
in p0,
in p1,
in p2,
in p3 )

Given the cubic polynomial defined with Bezier polygon return the Hermite data for the same polynomial

◆ cut()

function BaseHermite::cut ( in ignoredArg,
in varargin )

Cut the cubic polynomial \( \mathbf{p}(t) \) defined on Hermite data on the interval [a,b] and return the new Hermite data

[new_P0,new_P1,new_T0,new_T1] = object.cut( a, b, P0, P1, T0, T1 )
[new_P0,new_P1,new_T0,new_T1] = object.cut( a, b, P0, P1, T0, T1, H )

The parametrization of the new Hermite data is on [0,1]

◆ eval()

function BaseHermite::eval ( in ignoredArg,
in varargin )

Evaluate the cubic polynomial defined on Hermite data:

\begin{eqnarray*} \mathbf{p}(t) &=& h_0(t)\mathbf{p}_0+ h_1(t)\mathbf{p}_1+ h_2(t)\mathbf{t}_0+ h_3(t)\mathbf{t}_1 \\ \mathbf{p}(t,H) &=& h_0(t/H)\mathbf{p}_0+ h_1(t/H)\mathbf{p}_1+ H (h_2(t/H)\mathbf{t}_0+ h_3(t/H)\mathbf{t}_1) \end{eqnarray*}

values = object.eval( t, P0, P1, T0, T1 )
values = object.eval( t, P0, P1, T0, T1, H )

◆ eval5()

function BaseHermite::eval5 ( in ignoredArg,
in varargin )

Evaluate the quintic polynomial defined on Hermite data:

\begin{eqnarray*} \mathbf{p}(t) &=& h_0(t)\mathbf{p}_0+ h_1(t)\mathbf{p}_1+ h_2(t)\mathbf{t}_0+ h_3(t)\mathbf{t}_1+ h_4(t)\mathbf{j}_0+ h_5(t)\mathbf{j}_1 \\ \mathbf{p}(t,H) &=& h_0(t/H)\mathbf{p}_0+ h_1(t/H)\mathbf{p}_1+ H (h_2(t/H)\mathbf{t}_0+ h_3(t/H)\mathbf{t}_1)+ H^2 (h_4(t/H)\mathbf{t}_0+ h_5(t/H)\mathbf{j}_1) \end{eqnarray*}

values = object.eval5( t, P0, P1, T0, T1, J0, J1 )
values = object.eval5( t, P0, P1, T0, T1, J0, J1, H )

◆ eval5_D()

function BaseHermite::eval5_D ( in ignoredArg,
in varargin )

Evaluate the derivative \( \mathbf{p}'(t) \) of the quintic polynomial defined on Hermite data:

values = object.eval5_D( t, P0, P1, T0, T1, J0, J1 )
values = object.eval5_D( t, P0, P1, T0, T1, J0, J1, H )

◆ eval5_DD()

function BaseHermite::eval5_DD ( in ignoredArg,
in varargin )

Evaluate the second derivative \( \mathbf{p}''(t) \) of the quintic polynomial defined on Hermite data:

values = object.eval5_DD( t, P0, P1, T0, T1, J0, J1 )
values = object.eval5_DD( t, P0, P1, T0, T1, J0, J1, H )

◆ eval5_DDD()

function BaseHermite::eval5_DDD ( in ignoredArg,
in varargin )

Evaluate the third derivative \( \mathbf{p}'''(t) \) of the quintic polynomial defined on Hermite data:

values = object.eval5_DDD( t, P0, P1, T0, T1, J0, J1 )
values = object.eval5_DDD( t, P0, P1, T0, T1, J0, J1, H )

◆ eval5_DDDD()

function BaseHermite::eval5_DDDD ( in ignoredArg,
in varargin )

Evaluate the 4th derivative \( \mathbf{p}''''(t) \) of the quintic polynomial defined on Hermite data:

values = object.eval5_DDDD( t, P0, P1, T0, T1, J0, J1 )
values = object.eval5_DDDD( t, P0, P1, T0, T1, J0, J1, H )

◆ eval5_DDDDD()

function BaseHermite::eval5_DDDDD ( in ignoredArg,
in varargin )

Evaluate the 5th derivative \( \mathbf{p}'''''(t) \) of the quintic polynomial defined on Hermite data:

values = object.eval5_DDDDD( t, P0, P1, T0, T1, J0, J1 )
values = object.eval5_DDDDD( t, P0, P1, T0, T1, J0, J1, H )

◆ eval_D()

function BaseHermite::eval_D ( in ignoredArg,
in varargin )

Evaluate the derivative \( \mathbf{p}'(t) \) of the cubic polynomial defined on Hermite data:

values = object.eval_D( t, P0, P1, T0, T1 )
values = object.eval_D( t, P0, P1, T0, T1, H )

◆ eval_DD()

function BaseHermite::eval_DD ( in ignoredArg,
in varargin )

Evaluate the second derivative \( \mathbf{p}''(t) \) of the cubic polynomial defined on Hermite data:

values = object.eval_DD( t, P0, P1, T0, T1 )
values = object.eval_DD( t, P0, P1, T0, T1, H )

◆ eval_DDD()

function BaseHermite::eval_DDD ( in ignoredArg,
in varargin )

Evaluate the third derivative \( \mathbf{p}'''(t) \) of the cubic polynomial defined on Hermite data:

values = object.eval_DDD( t, P0, P1, T0, T1 )
values = object.eval_DDD( t, P0, P1, T0, T1, H )

◆ hermite_to_bezier()

function BaseHermite::hermite_to_bezier ( in ignoredArg,
in p0,
in p1,
in t0,
in t1 )

Given the cubic polynomial defined on Hermite data:

\begin{eqnarray*} \mathbf{p}(t) = h_0(t)\mathbf{p}_0 + h_1(t)\mathbf{p}_1 + h_2(t)\mathbf{t}_0 + h_3(t)\mathbf{t}_1 \end{eqnarray*}

return the Bezier polynomial of the same cubic


The documentation for this class was generated from the following file:
  • /Users/enrico/Ricerca/develop/PINS/pins-mechatronix/LibSources/submodules/Splines/toolbox/lib/BaseHermite.m