Rigid Body Dynamics Library
SegmentedQuinticBezierToolkit Class Reference

#include <SegmentedQuinticBezierToolkit.h>

Static Public Member Functions

static double scaleCurviness (double curviness)
 
static double calcU (double ax, const RigidBodyDynamics::Math::VectorNd &bezierPtsX, double tol, int maxIter)
 
static int calcIndex (double x, const RigidBodyDynamics::Math::MatrixNd &bezierPtsX)
 
static int calcIndex (double x, const std::vector< RigidBodyDynamics::Math::VectorNd > &bezierPtsX)
 
static double calcQuinticBezierCurveVal (double u, const RigidBodyDynamics::Math::VectorNd &pts)
 
static double calcQuinticBezierCurveDerivU (double u, const RigidBodyDynamics::Math::VectorNd &pts, int order)
 
static double calcQuinticBezierCurveDerivDYDX (double u, const RigidBodyDynamics::Math::VectorNd &xpts, const RigidBodyDynamics::Math::VectorNd &ypts, int order)
 
static RigidBodyDynamics::Math::MatrixNd calcQuinticBezierCornerControlPoints (double x0, double y0, double dydx0, double x1, double y1, double dydx1, double curviness)
 

Static Private Member Functions

static void printMatrixToFile (const RigidBodyDynamics::Math::VectorNd &col0, const RigidBodyDynamics::Math::MatrixNd &data, std::string &filename)
 
static void printBezierSplineFitCurves (const Function_< double > &curveFit, RigidBodyDynamics::Math::MatrixNd &ctrlPts, RigidBodyDynamics::Math::VectorNd &xVal, RigidBodyDynamics::Math::VectorNd &yVal, std::string &filename)
 
static double clampU (double u)
 

Detailed Description

Definition at line 177 of file SegmentedQuinticBezierToolkit.h.

Member Function Documentation

◆ calcIndex() [1/2]

static int calcIndex ( double  x,
const RigidBodyDynamics::Math::MatrixNd bezierPtsX 
)
static

Given a set of Bezier curve control points, return the index of the set of control points that x lies within.

Parameters
xA value that is interpolated by the set of Bezier curves
bezierPtsXA matrix of 6xn Bezier control points

aborts -If the index is not located within this set of Bezier points

Given a set of Bezier curve control points, return the index of the set of control points that x lies within. This function has been coded assuming a small number of Bezier curve sets (less than 10), and so, it simply scans through the Bezier curve sets until it finds the correct one.

Computational Costs Quoted for a Bezier curve set containing 1 to 5 curves.

  ~9-25

Example:

//The first set of spline points
mX(0,0) = -2;
mX(1,0) = -1;
mX(2,0) = -1;
mX(3,0) = 1;
mX(4,0) = 1;
mX(5,0) = 2;
//The second set of spline points
mX(0,1) = 2;
mX(1,1) = 3;
mX(2,1) = 3;
mX(3,1) = 5;
mX(4,1) = 5;
mX(5,1) = 6;
//The value of x for which we want the index for
double xVal = 1.75;
static int calcIndex(double x, const RigidBodyDynamics::Math::MatrixNd &bezierPtsX)

◆ calcIndex() [2/2]

static int calcIndex ( double  x,
const std::vector< RigidBodyDynamics::Math::VectorNd > &  bezierPtsX 
)
static

◆ calcQuinticBezierCornerControlPoints()

static RigidBodyDynamics::Math::MatrixNd calcQuinticBezierCornerControlPoints ( double  x0,
double  y0,
double  dydx0,
double  x1,
double  y1,
double  dydx1,
double  curviness 
)
static

Calculates the location of quintic Bezier curve control points to create a C shaped curve like that shown in the figure. Using a series of these simple and predictably shaped Bezier curves it is easy to build quite complex curves.

Parameters
x0First intercept x location
y0First intercept y location
dydx0First intercept slope
x1Second intercept x location
y1Second intercept y location
dydx1Second intercept slope
curvinessA parameter that ranges between 0 and 1 to denote a straight line or a curve

aborts -If the curviness parameter is less than 0, or greater than 1; -If the points and slopes are chosen so that an "S" shaped curve would be produced. This is tested by examining the points (x0,y0) and (x1,y1) together with the intersection (xC,yC) of the lines beginning at these points with slopes of dydx0 and dydx1 form a triangle. If the line segment from (x0,y0) to (x1,y1) is not the longest line segment, an exception is thrown. This is an overly conservative test as it prevents very deep 'V' shapes from being respresented.

Returns
a RigidBodyDynamics::Math::MatrixNd of 6 points Matrix(6,2) that correspond to the X, and Y control points for a quintic Bezier curve that has the above properties

Calculates the location of quintic Bezier curve control points to create a C shaped curve that intersects points 0 (x0, y0) and point 1 (x1, y1) with slopes dydx0 and dydx1 respectively, and a second derivative of 0. The curve that results can approximate a line (curviness = 0), or in a smooth C shaped curve (curviniess = 1)

The current implementation of this function is not optimized in anyway and has the following costs:

Computational Costs

  ~55 flops

Example:

double x0 = 1;
double y0 = 0;
double dydx0 = 0;
double x1 = 1.04;
double y1 = 1;
double dydx1 = 43;
double c = 0.75;
calcQuinticBezierCornerControlPoints(x0, y0, dydx0,x1,y1,dydx01,
c);
static RigidBodyDynamics::Math::MatrixNd calcQuinticBezierCornerControlPoints(double x0, double y0, double dydx0, double x1, double y1, double dydx1, double curviness)

◆ calcQuinticBezierCurveDerivDYDX()

static double calcQuinticBezierCurveDerivDYDX ( double  u,
const RigidBodyDynamics::Math::VectorNd xpts,
const RigidBodyDynamics::Math::VectorNd ypts,
int  order 
)
static

Calculates the value of dydx of a quintic Bezier curve derivative at u.

Parameters
uThe u value of interest. Note that u must be [0,1]
xptsThe 6 control points associated with the x axis
yptsThe 6 control points associated with the y axis
orderThe order of the derivative. Currently only orders from 1-6 can be evaluated

aborts -If u is outside [0,1] -If xpts is not 6 elements long -If ypts is not 6 elements long -If the order is less than 1 -If the order is greater than 6

Return values
Thevalue of (d^n y)/(dx^n) evaluated at u

Calculates the value of dydx of a quintic Bezier curve derivative at u. Note that a 2D Bezier curve can have an infinite number of derivatives, because x and y are functions of u. Thus

dy/dx = (dy/du)/(dx/du)
d^2y/dx^2 = d/du(dy/dx)*du/dx
    = [(d^2y/du^2)*(dx/du) - (dy/du)*(d^2x/du^2)]/(dx/du)^2 
    *(1/(dx/du))
etc.

Computational Costs

This obviously only functions when the Bezier curve in question has a finite derivative. Additionally, higher order derivatives are more numerically expensive to evaluate than lower order derivatives. For example, here are the number of operations required to compute the following derivatives

  Name  : flops
  dy/dx   : ~102
  d2y/dx2 : ~194
  d3y/dx3 : ~321
  d4y/dx4 : ~426
  d5y/dx5 : ~564
  d6y/dx6 : ~739

Example:

double u = 0.5;
vX(0) = 1;
vX(1) = 1.01164;
vX(2) = 1.01164;
vX(3) = 1.02364;
vX(4) = 1.02364;
vY(5) = 1.04;
vY(0) = 0;
vY(1) = 3e-16;
vY(2) = 3e-16;
vY(3) = 0.3;
vY(4) = 0.3;
vY(5) = 1;
u,vX, vY, 2);
static double calcQuinticBezierCurveDerivDYDX(double u, const RigidBodyDynamics::Math::VectorNd &xpts, const RigidBodyDynamics::Math::VectorNd &ypts, int order)

◆ calcQuinticBezierCurveDerivU()

static double calcQuinticBezierCurveDerivU ( double  u,
const RigidBodyDynamics::Math::VectorNd pts,
int  order 
)
static

Calculates the value of a quintic Bezier derivative curve at value u.

Parameters
uThe independent variable of a Bezier curve, which ranges between 0.0 and 1.0.
ptsThe locations of the control points in 1 dimension.
orderThe desired order of the derivative. Order must be >= 1

aborts -u is outside [0,1] -pts is not 6 elements long -if order is less than 1

Returns
The value of du/dx of Bezier curve located at u.

Calculates the value of a quintic Bezier derivative curve at value u. This calculation is acheived by taking the derivative of the row vector uV and multiplying it by the 6x6 coefficient matrix associated with a quintic Bezier curve, by the vector of Bezier control points, pV, in a particular dimension.

Pseudo code for the first derivative (order == 1) would be

  uV = [5*u^4 4*u^3 3*u^2 2u 1 0];

  cM = [ -1     5   -10    10  -5   1;
          5   -20    30   -20   5   0;
        -10    30   -30    10   0   0;
         10   -20    10     0   0   0;
         -5     5     0     0   0   0;
          1     0     0     0   0   0 ];
  pV = [x1; x2; x3; x4; x5; x6];

  dxdu = (uV*cM)*pV

Note that the derivative of uV only needed to be computed to compute dxdu. This process is continued for all 5 derivatives of x(u) until the sixth and all following derivatives, which are 0. Higher derivatives w.r.t. to U are less expensive to compute than lower derivatives.

Computational Costs

  dy/dx  : ~50 flops
  d2x/du2: ~43 flops
  d3x/du3: ~34 flops
  d4x/du4: ~26 flops
  d5x/du5: ~15 flops
  d6x/du6: ~1  flop

Example:

double u = 0.5;
//Choose the control points
vX(0) = -2;
vX(1) = 0;
vX(2) = 0;
vX(3) = 4;
vX(4) = 4;
vX(5) = 6;
double dxdu =calcQuinticBezierCurveDerivU(u,vX,1);
static double calcQuinticBezierCurveDerivU(double u, const RigidBodyDynamics::Math::VectorNd &pts, int order)

◆ calcQuinticBezierCurveVal()

static double calcQuinticBezierCurveVal ( double  u,
const RigidBodyDynamics::Math::VectorNd pts 
)
static

Calculates the value of a quintic Bezier curve at value u.

Parameters
uThe independent variable of a Bezier curve, which ranges between 0.0 and 1.0.
ptsThe locations of the control points in 1 dimension.

aborts -If u is outside of [0,1] -if pts has a length other than 6

Returns
The value of the Bezier curve located at u.

Calculates the value of a quintic Bezier curve at value u. This calculation is acheived by mulitplying a row vector comprised of powers of u, by the 6x6 coefficient matrix associated with a quintic Bezier curve, by the vector of Bezier control points, pV, in a particular dimension. The code to compute the value of a quintic bezier curve has been optimized to have the following cost:

Computational Costs

 ~54 flops

The math this function executes is decribed in pseudo code as the following:

  uV = [u^5 u^4 u^3 u^2 u 1];

  cM = [ -1   5   -10   10  -5   1; 
          5  -20   30  -20   5   0; 
        -10   30  -30   10   0   0; 
         10  -20   10    0   0   0; 
         -5    5    0    0   0   0;
          1    0    0    0   0   0 ];
  pV = [x1; x2; x3; x4; x5; x6];

  xB = (uV*cM)*pV

Example:

double u = 0.5;
//Choose the control points
vX(0) = -2;
vX(1) = 0;
vX(2) = 0;
vX(3) = 4;
vX(4) = 4;
vX(5) = 6;
static double calcQuinticBezierCurveVal(double u, const RigidBodyDynamics::Math::VectorNd &pts)

◆ calcU()

static double calcU ( double  ax,
const RigidBodyDynamics::Math::VectorNd bezierPtsX,
double  tol,
int  maxIter 
)
static

This function will compute the u value that correesponds to the given x for a quintic Bezier curve.

Parameters
axThe x value
bezierPtsXThe 6 Bezier point values
tolThe desired tolerance on u.
maxIterThe maximum number of Newton iterations allowed

aborts -if ax is outside the range defined in this Bezier spline section -if the desired tolerance is not met -if the derivative goes to 0 to machine precision

This function will compute the u value that correesponds to the given x for a quintic Bezier curve. This is accomplished by using an approximate spline inverse of u(x) to get a very good initial guess, and then one or two Newton iterations to polish the answer to the desired tolerance.

Computational Costs

  ~219 flops

Example:

double xVal = 2;
//Choose the control points
vX(0) = -2;
vX(1) = 0;
vX(2) = 0;
vX(3) = 4;
vX(4) = 4;
vX(5) = 6;
//Now evalutate u at the given xVal
calcU(xVal,vX, 1e-12,20);
static double calcU(double ax, const RigidBodyDynamics::Math::VectorNd &bezierPtsX, double tol, int maxIter)

◆ clampU()

static double clampU ( double  u)
staticprivate

This function will return a value that is equal to u, except when u is outside of[0,1], then it is clamped to be 0, or 1

Parameters
uThe parameter to be clamped
Return values
ubut restricted to 0,1.

◆ printBezierSplineFitCurves()

static void printBezierSplineFitCurves ( const Function_< double > &  curveFit,
RigidBodyDynamics::Math::MatrixNd ctrlPts,
RigidBodyDynamics::Math::VectorNd xVal,
RigidBodyDynamics::Math::VectorNd yVal,
std::string &  filename 
)
staticprivate
Parameters
curveFita function that evaluates a curve
ctrlPtscontrol point locations for the fitted Bezier curve
xValthe x values at the control points
yValthe y values at the control points
filenameof the output file.

◆ printMatrixToFile()

static void printMatrixToFile ( const RigidBodyDynamics::Math::VectorNd col0,
const RigidBodyDynamics::Math::MatrixNd data,
std::string &  filename 
)
staticprivate

This function will print cvs file of the column vector col0 and the matrix data.

Parameters
col0A vector that must have the same number of rows as the data matrix. This column vector is printed as the first column
dataA matrix of data
filenameThe name of the file to print

◆ scaleCurviness()

static double scaleCurviness ( double  curviness)
static

This scales the users value of curviness to be between [0+delta, 1-delta] because if curviness is allowed to equal 0 or 1, the second derivative becomes quite violent and the resulting curve is difficult to fit splines to.

Parameters
curviness
Return values
ascaled version of curviness

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