YARP
Yet Another Robot Platform
yarp::math Namespace Reference

Namespaces

 impl
 

Classes

class  FrameTransform
 
class  NormRand
 A static class grouping function for normal random number generator. More...
 
class  Quaternion
 
class  Rand
 A static class grouping function for uniform random number generator. More...
 
class  RandnScalar
 A random number generator, normal distribution. More...
 
class  RandScalar
 A random number generator, uniform in the range 0-1. More...
 
class  Vec2D
 

Functions

yarp::sig::Matrix pile (const yarp::sig::Matrix &m1, const yarp::sig::Matrix &m2)
 Matrix-Matrix concatenation by column (defined in Math.h). More...
 
yarp::sig::Matrix pile (const yarp::sig::Vector &v, const yarp::sig::Matrix &m)
 Vector-Matrix concatenation (defined in Math.h). More...
 
yarp::sig::Matrix pile (const yarp::sig::Matrix &m, const yarp::sig::Vector &v)
 Matrix-Vector concatenation (defined in Math.h). More...
 
yarp::sig::Matrix pile (const yarp::sig::Vector &v1, const yarp::sig::Vector &v2)
 Vector-Vector concatenation (defined in Math.h). More...
 
yarp::sig::Matrix cat (const yarp::sig::Matrix &m1, const yarp::sig::Matrix &m2)
 Matrix-Matrix concatenation by row (defined in Math.h). More...
 
yarp::sig::Matrix cat (const yarp::sig::Matrix &m, const yarp::sig::Vector &v)
 Matrix-Vector concatenation (defined in Math.h). More...
 
yarp::sig::Matrix cat (const yarp::sig::Vector &v, const yarp::sig::Matrix &m)
 Vector-Matrix concatenation (defined in Math.h). More...
 
yarp::sig::Vector cat (const yarp::sig::Vector &v1, const yarp::sig::Vector &v2)
 Vector-Vector concatenation (defined in Math.h). More...
 
yarp::sig::Vector cat (const yarp::sig::Vector &v, double s)
 Vector-scalar concatenation (defined in Math.h). More...
 
yarp::sig::Vector cat (double s, const yarp::sig::Vector &v)
 Scalar-vector concatenation (defined in Math.h). More...
 
yarp::sig::Vector cat (double s1, double s2)
 Scalar-scalar concatenation (defined in Math.h). More...
 
yarp::sig::Vector cat (double s1, double s2, double s3)
 
yarp::sig::Vector cat (double s1, double s2, double s3, double s4)
 
yarp::sig::Vector cat (double s1, double s2, double s3, double s4, double s5)
 
double dot (const yarp::sig::Vector &a, const yarp::sig::Vector &b)
 Scalar product between vectors (defined in Math.h). More...
 
yarp::sig::Matrix outerProduct (const yarp::sig::Vector &a, const yarp::sig::Vector &b)
 Outer product between vectors (defined in Math.h). More...
 
yarp::sig::Vector cross (const yarp::sig::Vector &a, const yarp::sig::Vector &b)
 Compute the cross product between two vectors (defined in Math.h). More...
 
yarp::sig::Matrix crossProductMatrix (const yarp::sig::Vector &v)
 Compute the cross product matrix, that is a 3-by-3 skew-symmetric matrix (defined in Math.h). More...
 
bool crossProductMatrix (const yarp::sig::Vector &v, yarp::sig::Matrix &res)
 Compute the cross product matrix, that is a 3-by-3 skew-symmetric matrix (defined in Math.h). More...
 
double norm (const yarp::sig::Vector &v)
 Returns the Euclidean norm of the vector (defined in Math.h). More...
 
double norm2 (const yarp::sig::Vector &v)
 Returns the Euclidean squared norm of the vector (defined in Math.h). More...
 
double findMax (const yarp::sig::Vector &v)
 Returns the maximum of the elements of a real vector (defined in Math.h). More...
 
double findMin (const yarp::sig::Vector &v)
 Returns the minimum of the elements of a real vector (defined in Math.h). More...
 
yarp::sig::Vector zeros (int s)
 Creates a vector of zeros (defined in Math.h). More...
 
yarp::sig::Vector ones (int s)
 Creates a vector of ones (defined in Math.h). More...
 
yarp::sig::Matrix eye (int r, int c)
 Build an identity matrix (defined in Math.h). More...
 
yarp::sig::Matrix eye (int n)
 Build a square identity matrix (defined in Math.h). More...
 
yarp::sig::Matrix zeros (int r, int c)
 Build a matrix of zeros (defined in Math.h). More...
 
double det (const yarp::sig::Matrix &in)
 Computes the determinant of a matrix (defined in Math.h). More...
 
yarp::sig::Matrix luinv (const yarp::sig::Matrix &in)
 Invert a square matrix using LU-decomposition (defined in Math.h). More...
 
bool eigenValues (const yarp::sig::Matrix &in, yarp::sig::Vector &real, yarp::sig::Vector &img)
 Computes eigenvalues of the n-by-n real nonsymmetric matrix (defined in Math.h). More...
 
double sign (const double &v)
 Invert a symmetric and positive definite matrix using Cholesky decomposition (defined in Math.h). More...
 
yarp::sig::Vector sign (const yarp::sig::Vector &v)
 Returns the sign vector of a real vector, that is a vector with 1 if the value is positive, -1 if negative, 0 if equal to zero (defined in Math.h). More...
 
yarp::sig::Vector dcm2axis (const yarp::sig::Matrix &R)
 Converts a dcm (direction cosine matrix) rotation matrix R to axis/angle representation (defined in Math.h). More...
 
yarp::sig::Matrix axis2dcm (const yarp::sig::Vector &v)
 Returns a dcm (direction cosine matrix) rotation matrix R from axis/angle representation (defined in Math.h). More...
 
yarp::sig::Vector dcm2euler (const yarp::sig::Matrix &R)
 Converts a dcm (direction cosine matrix) rotation matrix to euler angles (ZYZ) (defined in Math.h). More...
 
yarp::sig::Matrix euler2dcm (const yarp::sig::Vector &euler)
 Converts euler angles (ZYZ) vector in the corresponding dcm (direction cosine matrix) rotation matrix (defined in Math.h). More...
 
yarp::sig::Vector dcm2rpy (const yarp::sig::Matrix &R)
 Converts a dcm (direction cosine matrix) rotation matrix to roll-pitch-yaw angles (defined in Math.h). More...
 
yarp::sig::Matrix rpy2dcm (const yarp::sig::Vector &rpy)
 Converts roll-pitch-yaw angles in the corresponding dcm (direction cosine matrix) rotation matrix (defined in Math.h). More...
 
yarp::sig::Vector dcm2ypr (const yarp::sig::Matrix &R)
 Converts a dcm (direction cosine matrix) rotation matrix to yaw-roll-pitch angles (defined in Math.h). More...
 
yarp::sig::Matrix ypr2dcm (const yarp::sig::Vector &ypr)
 Converts yaw-pitch-roll angles in the corresponding dcm (direction cosine matrix) rotation matrix (defined in Math.h). More...
 
yarp::sig::Matrix SE3inv (const yarp::sig::Matrix &H)
 Returns the inverse of a 4 by 4 rototranslational matrix (defined in Math.h). More...
 
yarp::sig::Matrix adjoint (const yarp::sig::Matrix &H)
 Returns the adjoint matrix of a given roto-translational matrix (defined in Math.h). More...
 
yarp::sig::Matrix adjointInv (const yarp::sig::Matrix &H)
 Returns the inverse of the adjoint matrix of a given roto-translational matrix (defined in Math.h). More...
 
void SVD (const yarp::sig::Matrix &in, yarp::sig::Matrix &U, yarp::sig::Vector &S, yarp::sig::Matrix &V)
 Factorize the M-by-N matrix 'in' into the singular value decomposition in = U S V^T (defined in SVD.h). More...
 
void SVDMod (const yarp::sig::Matrix &in, yarp::sig::Matrix &U, yarp::sig::Vector &S, yarp::sig::Matrix &V)
 Perform SVD decomposition on a MxN matrix (for M >= N) (defined in SVD.h). More...
 
void SVDJacobi (const yarp::sig::Matrix &in, yarp::sig::Matrix &U, yarp::sig::Vector &S, yarp::sig::Matrix &V)
 Perform SVD decomposition on a matrix using the Jacobi method (defined in SVD.h). More...
 
yarp::sig::Matrix pinv (const yarp::sig::Matrix &in, double tol=0.0)
 Perform the moore-penrose pseudo-inverse of a matrix (defined in SVD.h). More...
 
void pinv (const yarp::sig::Matrix &in, yarp::sig::Matrix &out, double tol=0.0)
 Perform the moore-penrose pseudo-inverse of a matrix (defined in SVD.h). More...
 
yarp::sig::Matrix pinv (const yarp::sig::Matrix &in, yarp::sig::Vector &sv, double tol=0.0)
 Perform the moore-penrose pseudo-inverse of a matrix (defined in SVD.h). More...
 
void pinv (const yarp::sig::Matrix &in, yarp::sig::Matrix &out, yarp::sig::Vector &sv, double tol=0.0)
 Perform the moore-penrose pseudo-inverse of a matrix (defined in SVD.h). More...
 
yarp::sig::Matrix pinvDamped (const yarp::sig::Matrix &in, yarp::sig::Vector &sv, double damp)
 Perform the damped pseudo-inverse of a matrix (defined in SVD.h). More...
 
yarp::sig::Matrix pinvDamped (const yarp::sig::Matrix &in, double damp)
 Perform the damped pseudo-inverse of a matrix (defined in SVD.h). More...
 
void pinvDamped (const yarp::sig::Matrix &in, yarp::sig::Matrix &out, double damp)
 Perform the damped pseudo-inverse of a matrix (defined in SVD.h). More...
 
void pinvDamped (const yarp::sig::Matrix &in, yarp::sig::Matrix &out, yarp::sig::Vector &sv, double damp)
 Perform the damped pseudo-inverse of a matrix (defined in SVD.h). More...
 
yarp::sig::Matrix projectionMatrix (const yarp::sig::Matrix &A, double tol=0.0)
 Compute the projection matrix of A, that is defined as A times its pseudoinverse: A*pinv(A) (defined in SVD.h). More...
 
void projectionMatrix (const yarp::sig::Matrix &A, yarp::sig::Matrix &out, double tol=0.0)
 Compute the projection matrix of A, that is defined as A times its pseudoinverse: A*pinv(A) (defined in SVD.h). More...
 
yarp::sig::Matrix nullspaceProjection (const yarp::sig::Matrix &A, double tol=0.0)
 Compute the nullspace projection matrix of A, that is defined as the difference between the identity matrix and the pseudoinverse of A times A: (I - pinv(A)*A) (defined in SVD.h). More...
 
void nullspaceProjection (const yarp::sig::Matrix &A, yarp::sig::Matrix &out, double tol=0.0)
 Compute the nullspace projection matrix of A, that is defined as the difference between the identity matrix and the pseudoinverse of A times A: (I - pinv(A)*A) (defined in SVD.h). More...
 

Function Documentation

◆ adjoint()

Matrix yarp::math::adjoint ( const yarp::sig::Matrix H)

Returns the adjoint matrix of a given roto-translational matrix (defined in Math.h).

The adjoint is a (6x6) matrix: [R , S(r)*R; 0, R] where R is the rotational part of H, and r the translational part.

Parameters
His the 4 by 4 rototranslational matrix.
Returns
the adjoint matrix

Definition at line 904 of file math.cpp.

◆ adjointInv()

Matrix yarp::math::adjointInv ( const yarp::sig::Matrix H)

Returns the inverse of the adjoint matrix of a given roto-translational matrix (defined in Math.h).

The inverse of an adjoint is a (6x6) matrix: [R^T , -S(R^T*r)*R^T; 0 , R^T] where R is the rotational part of H, and r the translational part.

Parameters
His the 4 by 4 rototranslational matrix.
Returns
the inverse of the adjoint matrix

Definition at line 930 of file math.cpp.

◆ axis2dcm()

Matrix yarp::math::axis2dcm ( const yarp::sig::Vector v)

Returns a dcm (direction cosine matrix) rotation matrix R from axis/angle representation (defined in Math.h).

Parameters
vis the axis/angle vector.
Returns
4 by 4 homogeneous matrix with the rotation components in the top left 3 by 3 submatrix.

Definition at line 689 of file math.cpp.

◆ cat() [1/10]

Matrix yarp::math::cat ( const yarp::sig::Matrix m,
const yarp::sig::Vector v 
)

Matrix-Vector concatenation (defined in Math.h).

Add a column at the end of a matrix.

Parameters
ma matrix nXm
va vector n
Returns
a matrix nX(m+1)

Definition at line 363 of file math.cpp.

◆ cat() [2/10]

Matrix yarp::math::cat ( const yarp::sig::Matrix m1,
const yarp::sig::Matrix m2 
)

Matrix-Matrix concatenation by row (defined in Math.h).

Parameters
m1a matrix nXm
m2a matrix nXr
Returns
a matrix nX(m+r)

Definition at line 349 of file math.cpp.

◆ cat() [3/10]

Matrix yarp::math::cat ( const yarp::sig::Vector v,
const yarp::sig::Matrix m 
)

Vector-Matrix concatenation (defined in Math.h).

Add a column at the beginning of a matrix.

Parameters
ma matrix nXm
va vector n
Returns
a matrix nX(m+1)

Definition at line 376 of file math.cpp.

◆ cat() [4/10]

Vector yarp::math::cat ( const yarp::sig::Vector v,
double  s 
)

Vector-scalar concatenation (defined in Math.h).

Create a vector by putting the scalar at the end of the vector.

Parameters
van n-vector
sa scalar
Returns
a (n+1)-vector

Definition at line 401 of file math.cpp.

◆ cat() [5/10]

Vector yarp::math::cat ( const yarp::sig::Vector v1,
const yarp::sig::Vector v2 
)

Vector-Vector concatenation (defined in Math.h).

Create a vector by putting two vectors side by side.

Parameters
v1an n-vector
v2an m-vector
Returns
a (n+m)-vector

Definition at line 389 of file math.cpp.

◆ cat() [6/10]

Vector yarp::math::cat ( double  s,
const yarp::sig::Vector v 
)

Scalar-vector concatenation (defined in Math.h).

Create a vector by putting the scalar at the beginning of the vector.

Parameters
sa scalar
van n-vector
Returns
a (n+1)-vector

Definition at line 412 of file math.cpp.

◆ cat() [7/10]

Vector yarp::math::cat ( double  s1,
double  s2 
)

Scalar-scalar concatenation (defined in Math.h).

Create a vector containing the two specified scalar values.

Parameters
s1a scalar
s2a scalar
Returns
a 2-vector

Definition at line 423 of file math.cpp.

◆ cat() [8/10]

Vector yarp::math::cat ( double  s1,
double  s2,
double  s3 
)

Definition at line 431 of file math.cpp.

◆ cat() [9/10]

Vector yarp::math::cat ( double  s1,
double  s2,
double  s3,
double  s4 
)

Definition at line 440 of file math.cpp.

◆ cat() [10/10]

Vector yarp::math::cat ( double  s1,
double  s2,
double  s3,
double  s4,
double  s5 
)

Definition at line 450 of file math.cpp.

◆ cross()

Vector yarp::math::cross ( const yarp::sig::Vector a,
const yarp::sig::Vector b 
)

Compute the cross product between two vectors (defined in Math.h).

Parameters
afirst input vector
bsecond input vector
Returns
axb

Definition at line 480 of file math.cpp.

◆ crossProductMatrix() [1/2]

Matrix yarp::math::crossProductMatrix ( const yarp::sig::Vector v)

Compute the cross product matrix, that is a 3-by-3 skew-symmetric matrix (defined in Math.h).

Parameters
vthe vector
Returns
the cross product matrix

Definition at line 491 of file math.cpp.

◆ crossProductMatrix() [2/2]

bool yarp::math::crossProductMatrix ( const yarp::sig::Vector v,
yarp::sig::Matrix res 
)

Compute the cross product matrix, that is a 3-by-3 skew-symmetric matrix (defined in Math.h).

Parameters
vthe vector
resthe cross product matrix
Returns
true if operation succeeded, false otherwise

Definition at line 504 of file math.cpp.

◆ dcm2axis()

Vector yarp::math::dcm2axis ( const yarp::sig::Matrix R)

Converts a dcm (direction cosine matrix) rotation matrix R to axis/angle representation (defined in Math.h).

Parameters
Ris the input matrix.
Returns
4 by 1 vector for the axis/angle representation.

Definition at line 647 of file math.cpp.

◆ dcm2euler()

Vector yarp::math::dcm2euler ( const yarp::sig::Matrix R)

Converts a dcm (direction cosine matrix) rotation matrix to euler angles (ZYZ) (defined in Math.h).

Three angles are returned in a vector with the following format:

\[ \mathbf{v} = [\alpha, \beta, \gamma ]\]

such that the returned matrix satisfies the following:

\[ R = R_z(\alpha) R_y(\beta) R_z(\gamma) \]

Parameters
Ris the input ZYZ rotation matrix.
Returns
3 by 1 vector for the Euler angles representation.

Definition at line 726 of file math.cpp.

◆ dcm2rpy()

Vector yarp::math::dcm2rpy ( const yarp::sig::Matrix R)

Converts a dcm (direction cosine matrix) rotation matrix to roll-pitch-yaw angles (defined in Math.h).

Three angles are returned in a vector with the following format:

\[ \mathbf{v} = [\psi, \theta, \phi ]\]

such that the returned matrix satisfies the following:

\[ R = R_z(\phi) R_y(\theta) R_x(\psi) \]

Parameters
Ris the input ZYX rotation matrix.
Returns
3 by 1 vector for the roll pitch-yaw-angles representation.

Definition at line 780 of file math.cpp.

◆ dcm2ypr()

Vector yarp::math::dcm2ypr ( const yarp::sig::Matrix R)

Converts a dcm (direction cosine matrix) rotation matrix to yaw-roll-pitch angles (defined in Math.h).

Three angles are returned in a vector with the following format:

\[ \mathbf{v} = [\phi, \theta, \psi ]\]

such that the returned matrix satisfies the following:

\[ R = R_x(\psi) R_y(\theta) R_z(\phi) \]

Parameters
Ris the input XYZ rotation matrix.
Returns
3 by 1 vector for the yaw-pitch-roll angles representation.

Definition at line 834 of file math.cpp.

◆ det()

double yarp::math::det ( const yarp::sig::Matrix in)

Computes the determinant of a matrix (defined in Math.h).

Parameters
inthe matrix

Definition at line 583 of file math.cpp.

◆ dot()

double yarp::math::dot ( const yarp::sig::Vector a,
const yarp::sig::Vector b 
)

Scalar product between vectors (defined in Math.h).

Returns
a^T*b, where a and b are column vectors

Definition at line 461 of file math.cpp.

◆ eigenValues()

bool yarp::math::eigenValues ( const yarp::sig::Matrix in,
yarp::sig::Vector real,
yarp::sig::Vector img 
)

Computes eigenvalues of the n-by-n real nonsymmetric matrix (defined in Math.h).

Parameters
innonsymmetric n-by-n matrix
realthe real part of eigen values
imgthe imaginary part of eigen values
Returns
the real and imaginary part of the eigen values in separate matrices

Definition at line 600 of file math.cpp.

◆ euler2dcm()

Matrix yarp::math::euler2dcm ( const yarp::sig::Vector euler)

Converts euler angles (ZYZ) vector in the corresponding dcm (direction cosine matrix) rotation matrix (defined in Math.h).

The three euler angles are specified in a vector with the following structure:

\[ \mathbf{v} = [\alpha, \beta, \gamma ]\]

and the returned matrix is:

\[ R = R_z(\alpha) R_y(\beta) R_z(\gamma) \]

Parameters
euleris the input vector (alpha=z-rotation, beta=y-rotation, gamma=z-rotation).
Returns
4 by 4 homogeneous matrix representing the ZYZ rotation with the rotation components in the top left 3 by 3 submatrix.

Definition at line 764 of file math.cpp.

◆ eye() [1/2]

Matrix yarp::math::eye ( int  n)

Build a square identity matrix (defined in Math.h).

Parameters
nnumber of rows and columns
Returns
the new matrix

Definition at line 570 of file math.cpp.

◆ eye() [2/2]

Matrix yarp::math::eye ( int  r,
int  c 
)

Build an identity matrix (defined in Math.h).

Parameters
rnumber of rows
cnumber of columns
Returns
the new matrix

Definition at line 562 of file math.cpp.

◆ findMax()

double yarp::math::findMax ( const yarp::sig::Vector v)

Returns the maximum of the elements of a real vector (defined in Math.h).

Parameters
vis the input vector.
Returns
max(v).

Definition at line 530 of file math.cpp.

◆ findMin()

double yarp::math::findMin ( const yarp::sig::Vector v)

Returns the minimum of the elements of a real vector (defined in Math.h).

Parameters
vis the input vector.
Returns
min(v).

Definition at line 541 of file math.cpp.

◆ luinv()

Matrix yarp::math::luinv ( const yarp::sig::Matrix in)

Invert a square matrix using LU-decomposition (defined in Math.h).

Parameters
insquare matrix
Returns
the inverse of the matrix

Definition at line 588 of file math.cpp.

◆ norm()

double yarp::math::norm ( const yarp::sig::Vector v)

Returns the Euclidean norm of the vector (defined in Math.h).

Parameters
vis the input vector.
Returns
||v||.

Definition at line 520 of file math.cpp.

◆ norm2()

double yarp::math::norm2 ( const yarp::sig::Vector v)

Returns the Euclidean squared norm of the vector (defined in Math.h).

Parameters
vis the input vector.
Returns
||v||^2.

Definition at line 525 of file math.cpp.

◆ nullspaceProjection() [1/2]

Matrix yarp::math::nullspaceProjection ( const yarp::sig::Matrix A,
double  tol = 0.0 
)

Compute the nullspace projection matrix of A, that is defined as the difference between the identity matrix and the pseudoinverse of A times A: (I - pinv(A)*A) (defined in SVD.h).

Multiplying this null projection matrix times a vector projects the vector in the nullspace of A.

Parameters
Ainput matrix
tolsingular values less than tol are set to zero
Returns
The projection matrix associated with the nullspace of A

Definition at line 178 of file SVD.cpp.

◆ nullspaceProjection() [2/2]

void yarp::math::nullspaceProjection ( const yarp::sig::Matrix A,
yarp::sig::Matrix out,
double  tol = 0.0 
)

Compute the nullspace projection matrix of A, that is defined as the difference between the identity matrix and the pseudoinverse of A times A: (I - pinv(A)*A) (defined in SVD.h).

Multiplying this projection matrix times a vector projects the vector in the range of A.

Parameters
Ainput matrix
outthe projection matrix associated with the nullspace of A
tolsingular values less than tol are set to zero

Definition at line 185 of file SVD.cpp.

◆ ones()

Vector yarp::math::ones ( int  s)

Creates a vector of ones (defined in Math.h).

Parameters
sthe size of the new vector
Returns
a copy of the new vector

Definition at line 557 of file math.cpp.

◆ outerProduct()

Matrix yarp::math::outerProduct ( const yarp::sig::Vector a,
const yarp::sig::Vector b 
)

Outer product between vectors (defined in Math.h).

Returns
a*b^T, where a and b are column vectors

Definition at line 469 of file math.cpp.

◆ pile() [1/4]

Matrix yarp::math::pile ( const yarp::sig::Matrix m,
const yarp::sig::Vector v 
)

Matrix-Vector concatenation (defined in Math.h).

Add a row at the end of a matrix.

Parameters
ma matrix nXm
va vector m
Returns
a matrix (n+1)Xm

Definition at line 311 of file math.cpp.

◆ pile() [2/4]

Matrix yarp::math::pile ( const yarp::sig::Matrix m1,
const yarp::sig::Matrix m2 
)

Matrix-Matrix concatenation by column (defined in Math.h).

Parameters
m1a matrix nXm
m2a matrix rXm
Returns
a matrix (n+r)Xm

Definition at line 297 of file math.cpp.

◆ pile() [3/4]

Matrix yarp::math::pile ( const yarp::sig::Vector v,
const yarp::sig::Matrix m 
)

Vector-Matrix concatenation (defined in Math.h).

Add a row at the beginning of a matrix.

Parameters
ma matrix nXm
va vector m
Returns
a matrix (n+1)Xm

Definition at line 324 of file math.cpp.

◆ pile() [4/4]

Matrix yarp::math::pile ( const yarp::sig::Vector v1,
const yarp::sig::Vector v2 
)

Vector-Vector concatenation (defined in Math.h).

Create a two row matrix by stacking two vectors.

Parameters
v1an n-vector
v2an n-vector
Returns
a matrix 2Xn

Definition at line 337 of file math.cpp.

◆ pinv() [1/4]

Matrix yarp::math::pinv ( const yarp::sig::Matrix in,
double  tol = 0.0 
)

Perform the moore-penrose pseudo-inverse of a matrix (defined in SVD.h).

Parameters
ininput matrix
tolsingular values less than tol are set to zero
Returns
pseudo-inverse of the matrix 'in'

Definition at line 49 of file SVD.cpp.

◆ pinv() [2/4]

void yarp::math::pinv ( const yarp::sig::Matrix in,
yarp::sig::Matrix out,
double  tol = 0.0 
)

Perform the moore-penrose pseudo-inverse of a matrix (defined in SVD.h).

Parameters
ininput matrix
outpseudo-inverse of the matrix 'in'
tolsingular values less than tol are set to zero

Definition at line 64 of file SVD.cpp.

◆ pinv() [3/4]

void yarp::math::pinv ( const yarp::sig::Matrix in,
yarp::sig::Matrix out,
yarp::sig::Vector sv,
double  tol = 0.0 
)

Perform the moore-penrose pseudo-inverse of a matrix (defined in SVD.h).

Parameters
ininput matrix
outpseudo-inverse of the matrix 'in'
svvector containing the singular values of the input matrix
tolsingular values less than tol are set to zero

Definition at line 96 of file SVD.cpp.

◆ pinv() [4/4]

Matrix yarp::math::pinv ( const yarp::sig::Matrix in,
yarp::sig::Vector sv,
double  tol = 0.0 
)

Perform the moore-penrose pseudo-inverse of a matrix (defined in SVD.h).

Parameters
ininput matrix
svvector containing the singular values of the input matrix
tolsingular values less than tol are set to zero
Returns
pseudo-inverse of the matrix 'in'

Definition at line 79 of file SVD.cpp.

◆ pinvDamped() [1/4]

Matrix yarp::math::pinvDamped ( const yarp::sig::Matrix in,
double  damp 
)

Perform the damped pseudo-inverse of a matrix (defined in SVD.h).

Parameters
ininput matrix
dampdamping factor

Definition at line 120 of file SVD.cpp.

◆ pinvDamped() [2/4]

void yarp::math::pinvDamped ( const yarp::sig::Matrix in,
yarp::sig::Matrix out,
double  damp 
)

Perform the damped pseudo-inverse of a matrix (defined in SVD.h).

Parameters
ininput matrix
outdamped pseudo-inverse of the matrix 'in'
dampdamping factor

Definition at line 129 of file SVD.cpp.

◆ pinvDamped() [3/4]

void yarp::math::pinvDamped ( const yarp::sig::Matrix in,
yarp::sig::Matrix out,
yarp::sig::Vector sv,
double  damp 
)

Perform the damped pseudo-inverse of a matrix (defined in SVD.h).

Parameters
ininput matrix
outdamped pseudo-inverse of the matrix 'in'
svvector containing the singular values of the input matrix
dampdamping factor

Definition at line 136 of file SVD.cpp.

◆ pinvDamped() [4/4]

Matrix yarp::math::pinvDamped ( const yarp::sig::Matrix in,
yarp::sig::Vector sv,
double  damp 
)

Perform the damped pseudo-inverse of a matrix (defined in SVD.h).

Parameters
ininput matrix
svvector containing the singular values of the input matrix
dampdamping factor

Definition at line 113 of file SVD.cpp.

◆ projectionMatrix() [1/2]

Matrix yarp::math::projectionMatrix ( const yarp::sig::Matrix A,
double  tol = 0.0 
)

Compute the projection matrix of A, that is defined as A times its pseudoinverse: A*pinv(A) (defined in SVD.h).

Multiplying this projection matrix times a vector projects the vector in the range of A.

Parameters
Ainput matrix
tolsingular values less than tol are set to zero
Returns
The projection matrix associated with the range of A

Definition at line 153 of file SVD.cpp.

◆ projectionMatrix() [2/2]

void yarp::math::projectionMatrix ( const yarp::sig::Matrix A,
yarp::sig::Matrix out,
double  tol = 0.0 
)

Compute the projection matrix of A, that is defined as A times its pseudoinverse: A*pinv(A) (defined in SVD.h).

Multiplying this projection matrix times a vector projects the vector in the range of A.

Parameters
Ainput matrix
outthe projection matrix associated with the range of A
tolsingular values less than tol are set to zero

Definition at line 160 of file SVD.cpp.

◆ rpy2dcm()

Matrix yarp::math::rpy2dcm ( const yarp::sig::Vector rpy)

Converts roll-pitch-yaw angles in the corresponding dcm (direction cosine matrix) rotation matrix (defined in Math.h).

The three angles are specified in a vector with the following structure:

\[ \mathbf{v} = [\psi, \theta, \phi ]\]

and the returned matrix is:

\[ R = R_z(\phi) R_y(\theta) R_x(\psi) \]

Parameters
rpyis the input vector (\psi=roll x-rotation,\theta=pitch y-rotation, \phi=yaw z-rotation).
Returns
4 by 4 homogeneous matrix representing the ZYX rotation with the rotation components in the top left 3 by 3 submatrix.

Definition at line 818 of file math.cpp.

◆ SE3inv()

Matrix yarp::math::SE3inv ( const yarp::sig::Matrix H)

Returns the inverse of a 4 by 4 rototranslational matrix (defined in Math.h).

Parameters
His the 4 by 4 rototranslational matrix.
Returns
inverse of 4 by 4 rototranslational matrix.
Note
about 5 times faster than pinv()

Definition at line 883 of file math.cpp.

◆ sign() [1/2]

double yarp::math::sign ( const double &  v)

Invert a symmetric and positive definite matrix using Cholesky decomposition (defined in Math.h).

Parameters
insymmetric and positive definite matrix
Returns
the inverse of the matrix Returns the sign of a real number: 1 if positive, -1 if negative, 0 if equal to zero (defined in Math.h)
Parameters
vis a real number.
Returns
sign(v).

Definition at line 632 of file math.cpp.

◆ sign() [2/2]

Vector yarp::math::sign ( const yarp::sig::Vector v)

Returns the sign vector of a real vector, that is a vector with 1 if the value is positive, -1 if negative, 0 if equal to zero (defined in Math.h).

Parameters
vis the input vector.
Returns
sign(v).

Definition at line 638 of file math.cpp.

◆ SVD()

void yarp::math::SVD ( const yarp::sig::Matrix in,
yarp::sig::Matrix U,
yarp::sig::Vector S,
yarp::sig::Matrix V 
)

Factorize the M-by-N matrix 'in' into the singular value decomposition in = U S V^T (defined in SVD.h).

The diagonal elements of the singular value matrix S are stored in the vector S. The singular values are non-negative and form a non-increasing sequence from S_1 to S_N. The matrix V contains the elements of V in untransposed form. To form the product U S V^T it is necessary to take the transpose of V. Defining K as min(M, N) the the input matrices are:

Parameters
ininput M-by-N matrix to decompose
Uoutput M-by-K orthogonal matrix
Soutput K-dimensional vector containing the diagonal entries of the diagonal matrix S
Voutput N-by-K orthogonal matrix
Note
If U, S, V do not have the expected sizes they are resized automatically.
The routine computes the thin version of the SVD. Mathematically, the full SVD is defined with U and V as square orthogonal matrices and S as an M-by-N diagonal matrix.
This function uses the Jacobi SVD algorithm.

Definition at line 25 of file SVD.cpp.

◆ SVDJacobi()

void yarp::math::SVDJacobi ( const yarp::sig::Matrix in,
yarp::sig::Matrix U,
yarp::sig::Vector S,
yarp::sig::Matrix V 
)

Perform SVD decomposition on a matrix using the Jacobi method (defined in SVD.h).

The Jacobi method can compute singular values to higher relative accuracy than Golub-Reinsch algorithms.

Note
If U, S, V do not have the expected sizes they are resized automatically.
This function uses the Jacobi SVD algorithm.

Definition at line 35 of file SVD.cpp.

◆ SVDMod()

void yarp::math::SVDMod ( const yarp::sig::Matrix in,
yarp::sig::Matrix U,
yarp::sig::Vector S,
yarp::sig::Matrix V 
)

Perform SVD decomposition on a MxN matrix (for M >= N) (defined in SVD.h).

Note
If U, S, V do not have the expected sizes they are resized automatically.
This function uses the Jacobi SVD algorithm.

Definition at line 30 of file SVD.cpp.

◆ ypr2dcm()

Matrix yarp::math::ypr2dcm ( const yarp::sig::Vector ypr)

Converts yaw-pitch-roll angles in the corresponding dcm (direction cosine matrix) rotation matrix (defined in Math.h).

The three angles are specified in a vector with the following structure:

\[ \mathbf{v} = [\phi, \theta, \psi ]\]

and the returned matrix is:

\[ R = R_x(\psi) R_y(\theta) R_z(\phi) \]

Parameters
ypris the input vector (\phi=yaw z-rotation, \theta=pitch y-rotation, \psi=roll x-rotation).
Returns
4 by 4 homogeneous matrix representing the XYZ rotation with the rotation components in the top left 3 by 3 submatrix.

Definition at line 867 of file math.cpp.

◆ zeros() [1/2]

Matrix yarp::math::zeros ( int  r,
int  c 
)

Build a matrix of zeros (defined in Math.h).

Parameters
rnumber of rows
cnumber of columns

Definition at line 575 of file math.cpp.

◆ zeros() [2/2]

Vector yarp::math::zeros ( int  s)

Creates a vector of zeros (defined in Math.h).

Parameters
sthe size of the new vector
Returns
a copy of the new vector

Definition at line 552 of file math.cpp.