Drake
drake::math Namespace Reference

Namespaces

 internal
 
 quaternion_test
 
 test
 

Classes

struct  AutoDiffToGradientMatrix
 
struct  AutoDiffToValueMatrix
 
class  BarycentricMesh
 Represents a multi-linear function (from vector inputs to vector outputs) by interpolating between points on a mesh using (triangular) barycentric interpolation. More...
 
struct  GetSubMatrixGradientArray
 
struct  GetSubMatrixGradientSingleElement
 
struct  Gradient
 
struct  Gradient< Derived, Nq, 1 >
 
struct  GrayCodesMatrix
 GrayCodesMatrix::type returns an Eigen matrix of integers. More...
 
struct  GrayCodesMatrix< Eigen::Dynamic >
 
struct  MatGradMult
 
struct  MatGradMultMat
 
class  RollPitchYaw
 This class represents the orientation between two arbitrary frames A and D associated with a Space-fixed (extrinsic) X-Y-Z rotation by "roll-pitch-yaw" angles [r, p, y], which is equivalent to a Body-fixed (intrinsic) Z-Y-X rotation by "yaw-pitch-roll" angles [y, p, r]. More...
 
class  RotationMatrix
 This class represents a 3x3 rotation matrix between two arbitrary frames A and B and helps ensure users create valid rotation matrices. More...
 
class  Transform
 This class represents a rigid transform between two frames, which can be regarded in two ways. More...
 

Typedefs

template<typename Derived , int Nq>
using AutoDiffMatrixType = Eigen::Matrix< Eigen::AutoDiffScalar< Eigen::Matrix< typename Derived::Scalar, Nq, 1 >>, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime, 0, Derived::MaxRowsAtCompileTime, Derived::MaxColsAtCompileTime >
 The appropriate AutoDiffScalar gradient type given the value type and the number of derivatives at compile time. More...
 

Functions

template<typename Derived >
AutoDiffToValueMatrix< Derived >::type autoDiffToValueMatrix (const Eigen::MatrixBase< Derived > &auto_diff_matrix)
 
template<typename Derived >
std::enable_if< !std::is_same< typename Derived::Scalar, double >::value, Eigen::Matrix< typename Derived::Scalar::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime, 0, Derived::MaxRowsAtCompileTime, Derived::MaxColsAtCompileTime > >::type DiscardGradient (const Eigen::MatrixBase< Derived > &auto_diff_matrix)
 B = DiscardGradient(A) enables casting from a matrix of AutoDiffScalars to AutoDiffScalar::Scalar type, explicitly throwing away any gradient information. More...
 
template<typename Derived >
std::enable_if< std::is_same< typename Derived::Scalar, double >::value, const Eigen::MatrixBase< Derived > & >::type DiscardGradient (const Eigen::MatrixBase< Derived > &matrix)
 
template<typename _Scalar , int _Dim, int _Mode, int _Options>
std::enable_if< !std::is_same< _Scalar, double >::value, Eigen::Transform< typename _Scalar::Scalar, _Dim, _Mode, _Options > >::type DiscardGradient (const Eigen::Transform< _Scalar, _Dim, _Mode, _Options > &auto_diff_transform)
 
template<typename _Scalar , int _Dim, int _Mode, int _Options>
std::enable_if< std::is_same< _Scalar, double >::value, const Eigen::Transform< _Scalar, _Dim, _Mode, _Options > & >::type DiscardGradient (const Eigen::Transform< _Scalar, _Dim, _Mode, _Options > &transform)
 
template<typename Derived , typename DerivedAutoDiff >
void initializeAutoDiff (const Eigen::MatrixBase< Derived > &val, Eigen::MatrixBase< DerivedAutoDiff > &auto_diff_matrix, Eigen::DenseIndex num_derivatives=Eigen::Dynamic, Eigen::DenseIndex deriv_num_start=0)
 Initialize a single autodiff matrix given the corresponding value matrix. More...
 
template<int Nq = Eigen::Dynamic, typename Derived >
AutoDiffMatrixType< Derived, Nq > initializeAutoDiff (const Eigen::MatrixBase< Derived > &mat, Eigen::DenseIndex num_derivatives=-1, Eigen::DenseIndex deriv_num_start=0)
 Initialize a single autodiff matrix given the corresponding value matrix. More...
 
template<typename Derived >
void resizeDerivativesToMatchScalar (Eigen::MatrixBase< Derived > &mat, const typename Derived::Scalar &scalar)
 Resize derivatives vector of each element of a matrix to to match the size of the derivatives vector of a given scalar. More...
 
template<typename... Deriveds>
std::tuple< AutoDiffMatrixType< Deriveds, internal::totalSizeAtCompileTime< Deriveds... >)>... > initializeAutoDiffTuple (const Eigen::MatrixBase< Deriveds > &...args)
 Given a series of Eigen matrices, create a tuple of corresponding AutoDiff matrices with values equal to the input matrices and properly initialized derivative vectors. More...
 
template<typename Derived >
AutoDiffToGradientMatrix< Derived >::type autoDiffToGradientMatrix (const Eigen::MatrixBase< Derived > &auto_diff_matrix, int num_variables=Eigen::Dynamic)
 
template<typename Derived , typename DerivedGradient , typename DerivedAutoDiff >
void initializeAutoDiffGivenGradientMatrix (const Eigen::MatrixBase< Derived > &val, const Eigen::MatrixBase< DerivedGradient > &gradient, Eigen::MatrixBase< DerivedAutoDiff > &auto_diff_matrix)
 Initializes an autodiff matrix given a matrix of values and gradient matrix. More...
 
template<typename Derived , typename DerivedGradient >
AutoDiffMatrixType< Derived, DerivedGradient::ColsAtCompileTime > initializeAutoDiffGivenGradientMatrix (const Eigen::MatrixBase< Derived > &val, const Eigen::MatrixBase< DerivedGradient > &gradient)
 Creates and initializes an autodiff matrix given a matrix of values and gradient matrix. More...
 
template<typename DerivedGradient , typename DerivedAutoDiff >
void gradientMatrixToAutoDiff (const Eigen::MatrixBase< DerivedGradient > &gradient, Eigen::MatrixBase< DerivedAutoDiff > &auto_diff_matrix)
 
template<typename Derived >
std::enable_if< !std::is_same< typename Derived::Scalar, double >::value, Eigen::Matrix< typename Derived::Scalar::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime, 0, Derived::MaxRowsAtCompileTime, Derived::MaxColsAtCompileTime > >::type DiscardZeroGradient (const Eigen::MatrixBase< Derived > &auto_diff_matrix, const typename Eigen::NumTraits< typename Derived::Scalar::Scalar >::Real &precision=Eigen::NumTraits< typename Derived::Scalar::Scalar >::dummy_precision())
 B = DiscardZeroGradient(A, precision) enables casting from a matrix of AutoDiffScalars to AutoDiffScalar::Scalar type, but first checking that the gradient matrix is empty or zero. More...
 
template<typename Derived >
std::enable_if< std::is_same< typename Derived::Scalar, double >::value, const Eigen::MatrixBase< Derived > & >::type DiscardZeroGradient (const Eigen::MatrixBase< Derived > &matrix, double precision=0.)
 
template<typename _Scalar , int _Dim, int _Mode, int _Options>
std::enable_if< !std::is_same< _Scalar, double >::value, Eigen::Transform< typename _Scalar::Scalar, _Dim, _Mode, _Options > >::type DiscardZeroGradient (const Eigen::Transform< _Scalar, _Dim, _Mode, _Options > &auto_diff_transform, const typename Eigen::NumTraits< typename _Scalar::Scalar >::Real &precision=Eigen::NumTraits< typename _Scalar::Scalar >::dummy_precision())
 
template<typename _Scalar , int _Dim, int _Mode, int _Options>
std::enable_if< std::is_same< _Scalar, double >::value, const Eigen::Transform< _Scalar, _Dim, _Mode, _Options > & >::type DiscardZeroGradient (const Eigen::Transform< _Scalar, _Dim, _Mode, _Options > &transform, double precision=0.)
 
Eigen::MatrixXd ContinuousAlgebraicRiccatiEquation (const Eigen::Ref< const Eigen::MatrixXd > &A, const Eigen::Ref< const Eigen::MatrixXd > &B, const Eigen::Ref< const Eigen::MatrixXd > &Q, const Eigen::LLT< Eigen::MatrixXd > &R_cholesky)
 This is functionally the same as ContinuousAlgebraicRiccatiEquation(A, B, Q, R). More...
 
Eigen::MatrixXd ContinuousAlgebraicRiccatiEquation (const Eigen::Ref< const Eigen::MatrixXd > &A, const Eigen::Ref< const Eigen::MatrixXd > &B, const Eigen::Ref< const Eigen::MatrixXd > &Q, const Eigen::Ref< const Eigen::MatrixXd > &R)
 Computes the unique stabilizing solution S to the continuous-time algebraic Riccati equation: More...
 
template<typename v_Type , typename DtB_v_Type , typename w_AB_Type >
Vector3< typename v_Type::Scalar > ConvertTimeDerivativeToOtherFrame (const Eigen::MatrixBase< v_Type > &v_E, const Eigen::MatrixBase< DtB_v_Type > &DtB_v_E, const Eigen::MatrixBase< w_AB_Type > &w_AB_E)
 Given ᴮd/dt(v) (the time derivative in frame B of an arbitrary 3D vector v) and given ᴬωᴮ (frame B's angular velocity in another frame A), this method computes ᴬd/dt(v) (the time derivative in frame A of v) by: ᴬd/dt(v) = ᴮd/dt(v) + ᴬωᴮ x v. More...
 
template<typename Derived >
drake::Matrix3< typename Derived::Scalar > VectorToSkewSymmetric (const Eigen::MatrixBase< Derived > &p)
 
Eigen::MatrixXd DiscreteAlgebraicRiccatiEquation (const Eigen::Ref< const Eigen::MatrixXd > &A, const Eigen::Ref< const Eigen::MatrixXd > &B, const Eigen::Ref< const Eigen::MatrixXd > &Q, const Eigen::Ref< const Eigen::MatrixXd > &R)
 DiscreteAlgebraicRiccatiEquation function computes the unique stabilizing solution X to the discrete-time algebraic Riccati equation:

\[ A'XA - X - A'XB(B'XB+R)^{-1}B'XA + Q = 0 \]

. More...

 
template<typename Derived >
std::vector< Eigen::Triplet< typename Derived::Scalar > > SparseMatrixToTriplets (const Derived &matrix)
 For a sparse matrix, return a vector of triplets, such that we can reconstruct the matrix using setFromTriplet function. More...
 
template<typename Derived >
void SparseMatrixToRowColumnValueVectors (const Derived &matrix, std::vector< Eigen::Index > &row_indices, std::vector< Eigen::Index > &col_indices, std::vector< typename Derived::Scalar > &val)
 For a sparse matrix, return the row indices, the column indices, and value of the non-zero entries. More...
 
Eigen::Matrix3Xd UniformPtsOnSphereFibonacci (int num_points)
 Deterministically generates approximate evenly distributed points on a unit sphere. More...
 
template<typename Derived >
Eigen::Matrix< typename Derived::Scalar, 4, 1 > expmap2quat (const Eigen::MatrixBase< Derived > &v)
 
template<typename DerivedQ >
Eigen::Matrix< typename DerivedQ::Scalar, 3, 1 > quat2expmap (const Eigen::MatrixBase< DerivedQ > &q)
 
template<typename Derived1 , typename Derived2 >
Eigen::Matrix< typename Derived1::Scalar, 3, 1 > closestExpmap (const Eigen::MatrixBase< Derived1 > &expmap1, const Eigen::MatrixBase< Derived2 > &expmap2)
 
template<typename DerivedQ , typename DerivedE >
void quat2expmapSequence (const Eigen::MatrixBase< DerivedQ > &quat, const Eigen::MatrixBase< DerivedQ > &quat_dot, Eigen::MatrixBase< DerivedE > &expmap, Eigen::MatrixBase< DerivedE > &expmap_dot)
 
template<std::size_t Size>
std::array< int, Size > intRange (int start)
 
template<typename Derived >
Derived::PlainObject transposeGrad (const Eigen::MatrixBase< Derived > &dX, typename Derived::Index rows_X)
 
template<typename DerivedA , typename DerivedB , typename DerivedDA , typename DerivedDB >
MatGradMultMat< DerivedA, DerivedB, DerivedDA >::type matGradMultMat (const Eigen::MatrixBase< DerivedA > &A, const Eigen::MatrixBase< DerivedB > &B, const Eigen::MatrixBase< DerivedDA > &dA, const Eigen::MatrixBase< DerivedDB > &dB)
 
template<typename DerivedDA , typename DerivedB >
MatGradMult< DerivedDA, DerivedB >::type matGradMult (const Eigen::MatrixBase< DerivedDA > &dA, const Eigen::MatrixBase< DerivedB > &B)
 
template<typename Derived >
Eigen::Matrix< typename Derived::Scalar, Eigen::Dynamic, Eigen::Dynamic > getSubMatrixGradient (const Eigen::MatrixBase< Derived > &dM, const std::vector< int > &rows, const std::vector< int > &cols, typename Derived::Index M_rows, int q_start=0, typename Derived::Index q_subvector_size=-1)
 
template<int QSubvectorSize, typename Derived , std::size_t NRows, std::size_t NCols>
GetSubMatrixGradientArray< QSubvectorSize, Derived, NRows, NCols >::type getSubMatrixGradient (const Eigen::MatrixBase< Derived > &dM, const std::array< int, NRows > &rows, const std::array< int, NCols > &cols, typename Derived::Index M_rows, int q_start=0, typename Derived::Index q_subvector_size=QSubvectorSize)
 
template<int QSubvectorSize, typename Derived >
GetSubMatrixGradientSingleElement< QSubvectorSize, Derived >::type getSubMatrixGradient (const Eigen::MatrixBase< Derived > &dM, int row, int col, typename Derived::Index M_rows, typename Derived::Index q_start=0, typename Derived::Index q_subvector_size=QSubvectorSize)
 
template<typename DerivedA , typename DerivedB >
void setSubMatrixGradient (Eigen::MatrixBase< DerivedA > &dM, const Eigen::MatrixBase< DerivedB > &dM_submatrix, const std::vector< int > &rows, const std::vector< int > &cols, typename DerivedA::Index M_rows, typename DerivedA::Index q_start=0, typename DerivedA::Index q_subvector_size=-1)
 
template<int QSubvectorSize, typename DerivedA , typename DerivedB , std::size_t NRows, std::size_t NCols>
void setSubMatrixGradient (Eigen::MatrixBase< DerivedA > &dM, const Eigen::MatrixBase< DerivedB > &dM_submatrix, const std::array< int, NRows > &rows, const std::array< int, NCols > &cols, typename DerivedA::Index M_rows, typename DerivedA::Index q_start=0, typename DerivedA::Index q_subvector_size=QSubvectorSize)
 
template<int QSubvectorSize, typename DerivedDM , typename DerivedDMSub >
void setSubMatrixGradient (Eigen::MatrixBase< DerivedDM > &dM, const Eigen::MatrixBase< DerivedDMSub > &dM_submatrix, int row, int col, typename DerivedDM::Index M_rows, typename DerivedDM::Index q_start=0, typename DerivedDM::Index q_subvector_size=QSubvectorSize)
 
int GrayCodeToInteger (const Eigen::Ref< const Eigen::VectorXi > &gray_code)
 Converts the Gray code to an integer. More...
 
template<int NumDigits = Eigen::Dynamic>
GrayCodesMatrix< NumDigits >::type CalculateReflectedGrayCodes (int num_digits=NumDigits)
 Returns a matrix whose i'th row is the Gray code for integer i. More...
 
template<int MaxChunkSize = 10, class F , class Arg >
decltype(auto) jacobian (F &&f, Arg &&x)
 Computes a matrix of AutoDiffScalars from which both the value and the Jacobian of a function

\[ f:\mathbb{R}^{n\times m}\rightarrow\mathbb{R}^{p\times q} \]

(f: R^n*m -> R^p*q) can be extracted. More...

 
template<int MaxChunkSizeOuter = 10, int MaxChunkSizeInner = 10, class F , class Arg >
decltype(auto) hessian (F &&f, Arg &&x)
 Computes a matrix of AutoDiffScalars from which the value, Jacobian, and Hessian of a function

\[ f:\mathbb{R}^{n\times m}\rightarrow\mathbb{R}^{p\times q} \]

(f: R^n*m -> R^p*q) can be extracted. More...

 
template<typename Derived >
bool IsSymmetric (const Eigen::MatrixBase< Derived > &matrix)
 Determines if a matrix is symmetric. More...
 
template<typename Derived >
bool IsSymmetric (const Eigen::MatrixBase< Derived > &matrix, const typename Derived::Scalar &precision)
 Determines if a matrix is symmetric based on whether the difference between matrix(i, j) and matrix(j, i) is smaller than precision for all i, j. More...
 
template<typename Derived >
drake::MatrixX< typename Derived::Scalar > ToSymmetricMatrixFromLowerTriangularColumns (const Eigen::MatrixBase< Derived > &lower_triangular_columns)
 Given a column vector containing the stacked columns of the lower triangular part of a square matrix, returning a symmetric matrix whose lower triangular part is the same as the original matrix. More...
 
template<int rows, typename Derived >
Eigen::Matrix< typename Derived::Scalar, rows, rows > ToSymmetricMatrixFromLowerTriangularColumns (const Eigen::MatrixBase< Derived > &lower_triangular_columns)
 Given a column vector containing the stacked columns of the lower triangular part of a square matrix, returning a symmetric matrix whose lower triangular part is the same as the original matrix. More...
 
template<typename Derived >
void NormalizeVector (const Eigen::MatrixBase< Derived > &x, typename Derived::PlainObject &x_norm, typename drake::math::Gradient< Derived, Derived::RowsAtCompileTime, 1 >::type *dx_norm=nullptr, typename drake::math::Gradient< Derived, Derived::RowsAtCompileTime, 2 >::type *ddx_norm=nullptr)
 Computes the normalized vector, optionally with its gradient and second derivative. More...
 
template<class T >
Matrix3< T > ComputeBasisFromAxis (int axis_index, const Vector3< T > &axis_W)
 Creates a right-handed local basis from a given axis. More...
 
Eigen::MatrixXd DecomposePSDmatrixIntoXtransposeTimesX (const Eigen::Ref< const Eigen::MatrixXd > &Y, double zero_tol)
 For a symmetric positive semidefinite matrix Y, decompose it into XᵀX, where the number of rows in X equals to the rank of Y. More...
 
std::pair< Eigen::MatrixXd, Eigen::MatrixXd > DecomposePositiveQuadraticForm (const Eigen::Ref< const Eigen::MatrixXd > &Q, const Eigen::Ref< const Eigen::VectorXd > &b, double c, double tol=0)
 Rewrite a quadratic form xᵀQx + bᵀx + c to (Rx+d)ᵀ(Rx+d) where RᵀR = Q Rᵀd = b / 2 Notice that this decomposition is not unique. More...
 
template<typename Derived >
 DRAKE_DEPRECATED ("This code is deprecated. Use Eigen's Quaternion constructor ""Quaternion(w, x, y, z) -- not a home-brew 4-element vectors") Eigen
 
template<typename Scalar >
Eigen::Quaternion< Scalar > ClosestQuaternion (const Eigen::Quaternion< Scalar > &q0, const Eigen::Quaternion< Scalar > &q1)
 Returns a unit quaternion that represents the same orientation as q1, and has the "shortest" geodesic distance on the unit sphere to q0. More...
 
template<typename Derived >
Vector4< typename Derived::Scalar > quatConjugate (const Eigen::MatrixBase< Derived > &q)
 
template<typename Derived1 , typename Derived2 >
Vector4< typename Derived1::Scalar > quatProduct (const Eigen::MatrixBase< Derived1 > &q1, const Eigen::MatrixBase< Derived2 > &q2)
 
template<typename DerivedQ , typename DerivedV >
Vector3< typename DerivedV::Scalar > quatRotateVec (const Eigen::MatrixBase< DerivedQ > &q, const Eigen::MatrixBase< DerivedV > &v)
 
template<typename Derived1 , typename Derived2 >
Vector4< typename Derived1::Scalar > quatDiff (const Eigen::MatrixBase< Derived1 > &q1, const Eigen::MatrixBase< Derived2 > &q2)
 
template<typename Derived1 , typename Derived2 , typename DerivedU >
Derived1::Scalar quatDiffAxisInvar (const Eigen::MatrixBase< Derived1 > &q1, const Eigen::MatrixBase< Derived2 > &q2, const Eigen::MatrixBase< DerivedU > &u)
 
template<typename Derived >
Derived::Scalar quatNorm (const Eigen::MatrixBase< Derived > &q)
 
template<typename Derived1 , typename Derived2 , typename Scalar >
Vector4< Scalar > Slerp (const Eigen::MatrixBase< Derived1 > &q1, const Eigen::MatrixBase< Derived2 > &q2, const Scalar &interpolation_parameter)
 Q = Slerp(q1, q2, f) Spherical linear interpolation between two quaternions This function uses the implementation given in Algorithm 8 of [1]. More...
 
template<typename Derived >
Matrix3< typename Derived::Scalar > quat2rotmat (const Eigen::MatrixBase< Derived > &v)
 Computes the rotation matrix from quaternion representation. More...
 
template<typename T >
bool is_quaternion_in_canonical_form (const Eigen::Quaternion< T > &quat)
 This function tests whether a quaternion is in "canonical form" meaning that it tests whether the quaternion [w, x, y, z] has a non-negative w value. More...
 
template<typename T >
Eigen::Quaternion< T > QuaternionToCanonicalForm (const Eigen::Quaternion< T > &quat)
 This function returns a quaternion in its "canonical form" meaning that it returns a quaternion [w, x, y, z] with a non-negative w. More...
 
template<typename T >
bool AreQuaternionsEqualForOrientation (const Eigen::Quaternion< T > &quat1, const Eigen::Quaternion< T > &quat2, const T tolerance)
 This function tests whether two quaternions represent the same orientation. More...
 
template<typename T >
Vector4< T > CalculateQuaternionDtFromAngularVelocityExpressedInB (const Eigen::Quaternion< T > &quat_AB, const Vector3< T > &w_AB_B)
 This function calculates a quaternion's time-derivative from its quaternion and angular velocity. More...
 
template<typename T >
Vector3< T > CalculateAngularVelocityExpressedInBFromQuaternionDt (const Eigen::Quaternion< T > &quat_AB, const Vector4< T > &quatDt)
 This function calculates angular velocity from a quaternion and its time- derivative. More...
 
template<typename T >
CalculateQuaternionDtConstraintViolation (const Eigen::Quaternion< T > &quat, const Vector4< T > &quatDt)
 This function calculates how well a quaternion and its time-derivative satisfy the quaternion time-derivative constraint specified in [Kane, 1983] Section 1.13, equations 12-13, page 59. More...
 
template<typename T >
bool IsQuaternionValid (const Eigen::Quaternion< T > &quat, const double tolerance)
 This function tests if a quaternion satisfies the quaternion constraint specified in [Kane, 1983] Section 1.3, equation 4, page 12, i.e., a quaternion [w, x, y, z] must satisfy: w^2 + x^2 + y^2 + z^2 = 1. More...
 
template<typename T >
bool IsBothQuaternionAndQuaternionDtOK (const Eigen::Quaternion< T > &quat, const Vector4< T > &quatDt, const double tolerance)
 This function tests if a quaternion satisfies the time-derivative constraint specified in [Kane, 1983] Section 1.13, equation 13, page 59. More...
 
template<typename T >
bool IsQuaternionAndQuaternionDtEqualAngularVelocityExpressedInB (const Eigen::Quaternion< T > &quat, const Vector4< T > &quatDt, const Vector3< T > &w_B, const double tolerance)
 This function tests if a quaternion and a quaternions time-derivative can calculate and match an angular velocity to within a tolerance. More...
 
template<class Generator >
Eigen::AngleAxisd UniformlyRandomAngleAxis (Generator *generator)
 Generates a rotation (in the axis-angle representation) that rotates a point on the unit sphere to another point on the unit sphere with a uniform distribution over the sphere. More...
 
template<class Generator >
Eigen::Quaterniond UniformlyRandomQuaternion (Generator *generator)
 Generates a rotation (in the quaternion representation) that rotates a point on the unit sphere to another point on the unit sphere with a uniform distribution over the sphere. More...
 
template<class Generator >
RotationMatrix< doubleUniformlyRandomRotationMatrix (Generator *generator)
 Generates a rotation (in the rotation matrix representation) that rotates a point on the unit sphere to another point on the unit sphere with a uniform distribution over the sphere. More...
 
template<class Generator >
Eigen::Vector3d UniformlyRandomRPY (Generator *generator)
 Generates a rotation (in the roll-pitch-yaw representation) that rotates a point on the unit sphere to another point on the unit sphere with a uniform distribution over the sphere. More...
 
template<typename T >
 DRAKE_DEPRECATED ("This code is deprecated per issue #8323. ""Use constructor RollPitchYaw(quaternion).") Vector3< T > QuaternionToSpaceXYZ(const Eigen
 (Deprecated), use math::RollPitchYaw(quaternion). More...
 
template<typename Derived >
 DRAKE_DEPRECATED ("This code is deprecated per issue #8323. ""Use RollPitchYaw::ToQuaternion().") Quaternion< typename Derived
 (Deprecated), use math::RollPitchYaw(rpy).ToQuaternion(). More...
 
template<typename Derived >
drake::math::Gradient< Eigen::Matrix< typename Derived::Scalar, 3, 3 >, drake::kQuaternionSize >::type dquat2rotmat (const Eigen::MatrixBase< Derived > &quaternion)
 Computes the gradient of the function that converts a unit length quaternion to a rotation matrix. More...
 
template<typename DerivedR , typename DerivedDR >
drake::math::Gradient< Eigen::Matrix< typename DerivedR::Scalar, drake::kRpySize, 1 >, DerivedDR::ColsAtCompileTime >::type drotmat2rpy (const Eigen::MatrixBase< DerivedR > &R, const Eigen::MatrixBase< DerivedDR > &dR)
 Computes the gradient of the function that converts a rotation matrix to body-fixed z-y'-x'' Euler angles. More...
 
template<typename DerivedR , typename DerivedDR >
drake::math::Gradient< Eigen::Matrix< typename DerivedR::Scalar, drake::kQuaternionSize, 1 >, DerivedDR::ColsAtCompileTime >::type drotmat2quat (const Eigen::MatrixBase< DerivedR > &R, const Eigen::MatrixBase< DerivedDR > &dR)
 Computes the gradient of the function that converts rotation matrix to quaternion. More...
 
template<typename Derived >
Matrix3< typename Derived::Scalar > rpy2rotmat (const Eigen::MatrixBase< Derived > &rpy)
 
template<typename Derived >
 DRAKE_DEPRECATED ("This code is deprecated per issue #8323. ""Use RotationMatrix::ToQuaternion.") Vector4< typename Derived
 (Deprecated), use math::RotationMatrix::ToQuaternion More...
 
template<typename Derived >
Vector3< typename Derived::Scalar > rotmat2rpy (const Eigen::MatrixBase< Derived > &R)
 Computes SpaceXYZ Euler angles from rotation matrix. More...
 
template<typename Derived >
 DRAKE_DEPRECATED ("This code is deprecated per issue #8323. ""Use math::RotationMatrix(R) or math::rotmat2rpy(R) ""or math::RotationMatrix::ToQuaternionAsVector4.") VectorX< typename Derived
 (Deprecated), use math::RotationMatrix or math::rotmat2rpy or math::RotationMatrix::ToQuaternionAsVector4. More...
 
template<typename T >
Matrix3< T > XRotation (const T &theta)
 (Deprecated), use math::RotationMatrix::MakeXRotation More...
 
template<typename T >
Matrix3< T > YRotation (const T &theta)
 (Deprecated), use math::RotationMatrix::MakeYRotation More...
 
template<typename T >
Matrix3< T > ZRotation (const T &theta)
 (Deprecated), use math::RotationMatrix::MakeZRotation More...
 
template<typename Derived >
Matrix3< typename Derived::Scalar > ProjectMatToOrthonormalMat (const Eigen::MatrixBase< Derived > &M)
 (Deprecated), use math::RotationMatrix::ProjectToRotationMatrix More...
 
template<typename Derived >
Matrix3< typename Derived::Scalar > ProjectMatToRotMat (const Eigen::MatrixBase< Derived > &M)
 (Deprecated), use math::RotationMatrix::ProjectToRotationMatrix More...
 
template<typename Derived >
double ProjectMatToRotMatWithAxis (const Eigen::MatrixBase< Derived > &M, const Eigen::Ref< const Eigen::Vector3d > &axis, double angle_lb, double angle_ub)
 Projects a 3 x 3 matrix M onto SO(3). More...
 
template<class T1 , class T2 , class T3 >
T1 saturate (const T1 &value, const T2 &low, const T3 &high)
 Saturates the input value between upper and lower bounds. More...
 
template<class T1 , class T2 >
T1 wrap_to (const T1 &value, const T2 &low, const T2 &high)
 For variables that are meant to be periodic, (e.g. More...
 

Typedef Documentation

using AutoDiffMatrixType = Eigen::Matrix< Eigen::AutoDiffScalar<Eigen::Matrix<typename Derived::Scalar, Nq, 1>>, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime, 0, Derived::MaxRowsAtCompileTime, Derived::MaxColsAtCompileTime>

The appropriate AutoDiffScalar gradient type given the value type and the number of derivatives at compile time.

Function Documentation

bool drake::math::AreQuaternionsEqualForOrientation ( const Eigen::Quaternion< T > &  quat1,
const Eigen::Quaternion< T > &  quat2,
const T  tolerance 
)

This function tests whether two quaternions represent the same orientation.

This function converts each quaternion to its canonical form and tests whether the absolute value of the difference in corresponding elements of these canonical quaternions is within tolerance.

Parameters
quat1Quaternion [w, x, y, z] that relates two right-handed orthogonal unitary bases e.g., Ax, Ay, Az (A) to Bx, By, Bz (B). Note: quat is analogous to the rotation matrix R_AB.
quat2Quaternion with a description analogous to quat1.
toleranceNonnegative real scalar defining the allowable difference in the orientation described by quat1 and quat2.
Returns
true if quat1 and quat2 represent the same orientation (to within tolerance), otherwise false.

Here is the call graph for this function:

Here is the caller graph for this function:

AutoDiffToGradientMatrix<Derived>::type drake::math::autoDiffToGradientMatrix ( const Eigen::MatrixBase< Derived > &  auto_diff_matrix,
int  num_variables = Eigen::Dynamic 
)

Here is the call graph for this function:

Here is the caller graph for this function:

AutoDiffToValueMatrix<Derived>::type drake::math::autoDiffToValueMatrix ( const Eigen::MatrixBase< Derived > &  auto_diff_matrix)

Here is the caller graph for this function:

Vector3<T> drake::math::CalculateAngularVelocityExpressedInBFromQuaternionDt ( const Eigen::Quaternion< T > &  quat_AB,
const Vector4< T > &  quatDt 
)

This function calculates angular velocity from a quaternion and its time- derivative.

Algorithm from [Kane, 1983] Section 1.13, Pages 58-59.

Parameters
quat_ABQuaternion [w, x, y, z] that relates two right-handed orthogonal unitary bases e.g., Ax, Ay, Az (A) to Bx, By, Bz (B). Note: quat_AB is analogous to the rotation matrix R_AB.
quatDtTime-derivative of quat_AB, i.e. [ẇ, ẋ, ẏ, ż].
Return values
w_AB_BB's angular velocity in A, expressed in B.

Here is the caller graph for this function:

T drake::math::CalculateQuaternionDtConstraintViolation ( const Eigen::Quaternion< T > &  quat,
const Vector4< T > &  quatDt 
)

This function calculates how well a quaternion and its time-derivative satisfy the quaternion time-derivative constraint specified in [Kane, 1983] Section 1.13, equations 12-13, page 59.

For a quaternion [w, x, y, z], the quaternion must satisfy: w^2 + x^2 + y^2 + z^2 = 1, hence its time-derivative must satisfy: 2*(w*ẇ + x*ẋ + y*ẏ + z*ż) = 0.

Parameters
quatQuaternion [w, x, y, z] that relates two right-handed orthogonal unitary bases e.g., Ax, Ay, Az (A) to Bx, By, Bz (B). Note: A quaternion like quat_AB is analogous to the rotation matrix R_AB.
quatDtTime-derivative of quat, i.e., [ẇ, ẋ, ẏ, ż].
Return values
quaternionDt_constraint_violationThe amount the time- derivative of the quaternion constraint has been violated, which may be positive or negative (0 means the constraint is perfectly satisfied).

Here is the caller graph for this function:

Vector4<T> drake::math::CalculateQuaternionDtFromAngularVelocityExpressedInB ( const Eigen::Quaternion< T > &  quat_AB,
const Vector3< T > &  w_AB_B 
)

This function calculates a quaternion's time-derivative from its quaternion and angular velocity.

Algorithm from [Kane, 1983] Section 1.13, Pages 58-59.

Parameters
quat_ABQuaternion [w, x, y, z] that relates two right-handed orthogonal unitary bases e.g., Ax, Ay, Az (A) to Bx, By, Bz (B). Note: quat_AB is analogous to the rotation matrix R_AB.
w_AB_BB's angular velocity in A, expressed in B.
Return values
quatDtTime-derivative of quat_AB, i.e., [ẇ, ẋ, ẏ, ż].

Here is the caller graph for this function:

GrayCodesMatrix<NumDigits>::type drake::math::CalculateReflectedGrayCodes ( int  num_digits = NumDigits)

Returns a matrix whose i'th row is the Gray code for integer i.

Template Parameters
NumDigitsThe number of digits in the Gray code.
Parameters
num_digitsThe number of digits in the Gray code.
Returns
M. M is a matrix of size 2ᵏ x k, where k is num_digits. M.row(i) is the Gray code for integer i.

Here is the call graph for this function:

Here is the caller graph for this function:

Eigen::Matrix<typename Derived1::Scalar, 3, 1> drake::math::closestExpmap ( const Eigen::MatrixBase< Derived1 > &  expmap1,
const Eigen::MatrixBase< Derived2 > &  expmap2 
)

Here is the call graph for this function:

Here is the caller graph for this function:

Eigen::Quaternion<Scalar> drake::math::ClosestQuaternion ( const Eigen::Quaternion< Scalar > &  q0,
const Eigen::Quaternion< Scalar > &  q1 
)

Returns a unit quaternion that represents the same orientation as q1, and has the "shortest" geodesic distance on the unit sphere to q0.

Here is the caller graph for this function:

Matrix3<T> drake::math::ComputeBasisFromAxis ( int  axis_index,
const Vector3< T > &  axis_W 
)

Creates a right-handed local basis from a given axis.

Defines two other arbitrary axes such that the basis is orthonormal. The basis is R_WL, where W is the frame in which the input axis is expressed and L is a local basis such that v_W = R_WL * v_L.

Parameters
[in]axis_indexThe index of the axis (in the range [0,2]), with 0 corresponding to the x-axis, 1 corresponding to the y-axis, and z-corresponding to the z-axis.
[in]axis_WThe vector defining the basis's given axis expressed in frame W. The vector need not be a unit vector: this routine will normalize it.
Return values
R_WLThe computed basis.
Exceptions
std::logic_errorif the norm of axis_W is within 1e-10 to zero or axis_index does not lie in the range [0,2].

Here is the caller graph for this function:

Eigen::MatrixXd ContinuousAlgebraicRiccatiEquation ( const Eigen::Ref< const Eigen::MatrixXd > &  A,
const Eigen::Ref< const Eigen::MatrixXd > &  B,
const Eigen::Ref< const Eigen::MatrixXd > &  Q,
const Eigen::LLT< Eigen::MatrixXd > &  R_cholesky 
)

This is functionally the same as ContinuousAlgebraicRiccatiEquation(A, B, Q, R).

The Cholesky decomposition of R is passed in instead of R.

Here is the call graph for this function:

Here is the caller graph for this function:

Eigen::MatrixXd ContinuousAlgebraicRiccatiEquation ( const Eigen::Ref< const Eigen::MatrixXd > &  A,
const Eigen::Ref< const Eigen::MatrixXd > &  B,
const Eigen::Ref< const Eigen::MatrixXd > &  Q,
const Eigen::Ref< const Eigen::MatrixXd > &  R 
)

Computes the unique stabilizing solution S to the continuous-time algebraic Riccati equation:

///  S A + A' S - S B inv(R) B' S + Q = 0
/// 
@throws std::runtime_error if R is not positive definite.

Based on the Matrix Sign Function method outlined in this paper:
http://www.engr.iupui.edu/~skoskie/ECE684/Riccati_algorithms.pdf

Here is the call graph for this function:

Vector3<typename v_Type::Scalar> drake::math::ConvertTimeDerivativeToOtherFrame ( const Eigen::MatrixBase< v_Type > &  v_E,
const Eigen::MatrixBase< DtB_v_Type > &  DtB_v_E,
const Eigen::MatrixBase< w_AB_Type > &  w_AB_E 
)

Given ᴮd/dt(v) (the time derivative in frame B of an arbitrary 3D vector v) and given ᴬωᴮ (frame B's angular velocity in another frame A), this method computes ᴬd/dt(v) (the time derivative in frame A of v) by: ᴬd/dt(v) = ᴮd/dt(v) + ᴬωᴮ x v.

This mathematical operation is known as the "Transport Theorem" or the "Golden Rule for Vector Differentiation" [Mitiguy 2016, §7.3]. It was discovered by Euler in 1758. Its explicit notation with superscript frames was invented by Thomas Kane in 1950. Its use as the defining property of angular velocity was by Mitiguy in 1993.

In source code and comments, we use the following monogram notations: DtA_v = ᴬd/dt(v) denotes the time derivative in frame A of the vector v. DtA_v_E = [ᴬd/dt(v)]_E denotes the time derivative in frame A of vector v, with the resulting new vector quantity expressed in a frame E.

In source code, this mathematical operation is performed with all vectors expressed in the same frame E as [ᴬd/dt(v)]ₑ = [ᴮd/dt(v)]ₑ + [ᴬωᴮ]ₑ x [v]ₑ which in monogram notation is:

  DtA_v_E = DtB_v_E + w_AB_E x v_E

[Mitiguy 2016] Mitiguy, P., 2016. Advanced Dynamics & Motion Simulation.

std::pair< Eigen::MatrixXd, Eigen::MatrixXd > DecomposePositiveQuadraticForm ( const Eigen::Ref< const Eigen::MatrixXd > &  Q,
const Eigen::Ref< const Eigen::VectorXd > &  b,
double  c,
double  tol = 0 
)

Rewrite a quadratic form xᵀQx + bᵀx + c to (Rx+d)ᵀ(Rx+d) where RᵀR = Q Rᵀd = b / 2 Notice that this decomposition is not unique.

For example, with any permutation matrix P, we can define R₁ = P*R d₁ = P*d Then (R₁*x+d₁)ᵀ(R₁*x+d₁) gives the same quadratic form.

Parameters
QThe square matrix.
bThe vector containing the linear coefficients.
cThe constatnt term.
tolWe will determine if this quadratic form is always non-negative, by checking the Eigen values of the matrix [Q b/2] [bᵀ/2 c] are all greater than -tol.
Default: is 0.
Return values
(R,d).R and d have the same number of rows. R.cols() == x.rows(). The matrix X = [R d] has the same number of rows as the rank of
   [Q    b/2]
   [bᵀ/2   c]
Precondition
1. The quadratic form is always non-negative, namely the matrix
        [Q    b/2]
        [bᵀ/2   c]
        
is positive semidefinite.
  1. Q and b are of the correct size.
  2. tol is non-negative.
Exceptions
aruntime_error if the precondition is not satisfied.

Here is the call graph for this function:

Here is the caller graph for this function:

Eigen::MatrixXd DecomposePSDmatrixIntoXtransposeTimesX ( const Eigen::Ref< const Eigen::MatrixXd > &  Y,
double  zero_tol 
)

For a symmetric positive semidefinite matrix Y, decompose it into XᵀX, where the number of rows in X equals to the rank of Y.

Notice that this decomposition is not unique. For any orthonormal matrix U, s.t UᵀU = identity, X_prime = UX also satisfies X_primeᵀX_prime = Y. Here we only return one valid decomposition.

Parameters
YA symmetric positive semidefinite matrix.
zero_tolWe will need to check if some value (for example, the absolute value of Y's eigenvalues) is smaller than zero_tol. If it is, then we deem that value as 0.
Return values
X.The matrix X satisfies XᵀX = Y and X.rows() = rank(Y).
Precondition
1. Y is positive semidefinite.
  1. zero_tol is non-negative.
Exceptions
std::runtime_errorwhen the pre-conditions are not satisfied.

Here is the call graph for this function:

Here is the caller graph for this function:

std::enable_if< !std::is_same<typename Derived::Scalar, double>::value, Eigen::Matrix<typename Derived::Scalar::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime, 0, Derived::MaxRowsAtCompileTime, Derived::MaxColsAtCompileTime> >::type drake::math::DiscardGradient ( const Eigen::MatrixBase< Derived > &  auto_diff_matrix)

B = DiscardGradient(A) enables casting from a matrix of AutoDiffScalars to AutoDiffScalar::Scalar type, explicitly throwing away any gradient information.

For a matrix of type, e.g. MatrixX<AutoDiffXd> A, the comparable operation B = A.cast<double>() should (and does) fail to compile. Use DiscardGradient(A) if you want to force the cast (and explicitly declare that information is lost).

This method is overloaded to permit the user to call it for double types and AutoDiffScalar types (to avoid the calling function having to handle the two cases differently).

See also
DiscardZeroGradient

Here is the call graph for this function:

Here is the caller graph for this function:

std::enable_if< std::is_same<typename Derived::Scalar, double>::value, const Eigen::MatrixBase<Derived>&>::type drake::math::DiscardGradient ( const Eigen::MatrixBase< Derived > &  matrix)
std::enable_if< !std::is_same<_Scalar, double>::value, Eigen::Transform<typename _Scalar::Scalar, _Dim, _Mode, _Options> >::type drake::math::DiscardGradient ( const Eigen::Transform< _Scalar, _Dim, _Mode, _Options > &  auto_diff_transform)
See also
DiscardGradient().

Here is the call graph for this function:

std::enable_if<std::is_same<_Scalar, double>::value, const Eigen::Transform<_Scalar, _Dim, _Mode, _Options>&>::type drake::math::DiscardGradient ( const Eigen::Transform< _Scalar, _Dim, _Mode, _Options > &  transform)
See also
DiscardGradient().

Here is the call graph for this function:

std::enable_if< !std::is_same<typename Derived::Scalar, double>::value, Eigen::Matrix<typename Derived::Scalar::Scalar, Derived::RowsAtCompileTime, Derived::ColsAtCompileTime, 0, Derived::MaxRowsAtCompileTime, Derived::MaxColsAtCompileTime> >::type drake::math::DiscardZeroGradient ( const Eigen::MatrixBase< Derived > &  auto_diff_matrix,
const typename Eigen::NumTraits< typename Derived::Scalar::Scalar >::Real &  precision = Eigen::NumTraits<typename Derived::Scalar::Scalar>::dummy_precision() 
)

B = DiscardZeroGradient(A, precision) enables casting from a matrix of AutoDiffScalars to AutoDiffScalar::Scalar type, but first checking that the gradient matrix is empty or zero.

For a matrix of type, e.g. MatrixX<AutoDiffXd> A, the comparable operation B = A.cast<double>() should (and does) fail to compile. Use DiscardZeroGradient(A) if you want to force the cast (and the check).

This method is overloaded to permit the user to call it for double types and AutoDiffScalar types (to avoid the calling function having to handle the two cases differently).

Parameters
precisionis passed to Eigen's isZero(precision) to evaluate whether the gradients are zero.
Exceptions
std::runtime_errorif the gradients were not empty nor zero.
See also
DiscardGradient

Here is the call graph for this function:

Here is the caller graph for this function:

std::enable_if<std::is_same<typename Derived::Scalar, double>::value, const Eigen::MatrixBase<Derived>&>::type drake::math::DiscardZeroGradient ( const Eigen::MatrixBase< Derived > &  matrix,
double  precision = 0. 
)
See also
DiscardZeroGradient().

Here is the call graph for this function:

std::enable_if< !std::is_same<_Scalar, double>::value, Eigen::Transform<typename _Scalar::Scalar, _Dim, _Mode, _Options> >::type drake::math::DiscardZeroGradient ( const Eigen::Transform< _Scalar, _Dim, _Mode, _Options > &  auto_diff_transform,
const typename Eigen::NumTraits< typename _Scalar::Scalar >::Real &  precision = Eigen::NumTraits<typename _Scalar::Scalar>::dummy_precision() 
)
See also
DiscardZeroGradient().

Here is the call graph for this function:

std::enable_if< std::is_same<_Scalar, double>::value, const Eigen::Transform<_Scalar, _Dim, _Mode, _Options>&>::type drake::math::DiscardZeroGradient ( const Eigen::Transform< _Scalar, _Dim, _Mode, _Options > &  transform,
double  precision = 0. 
)
See also
DiscardZeroGradient().

Here is the call graph for this function:

Eigen::MatrixXd DiscreteAlgebraicRiccatiEquation ( const Eigen::Ref< const Eigen::MatrixXd > &  A,
const Eigen::Ref< const Eigen::MatrixXd > &  B,
const Eigen::Ref< const Eigen::MatrixXd > &  Q,
const Eigen::Ref< const Eigen::MatrixXd > &  R 
)

DiscreteAlgebraicRiccatiEquation function computes the unique stabilizing solution X to the discrete-time algebraic Riccati equation:

\[ A'XA - X - A'XB(B'XB+R)^{-1}B'XA + Q = 0 \]

.

Computes the unique stabilizing solution X to the discrete-time algebraic Riccati equation:

Exceptions
std::runtime_errorif Q is not positive semi-definite.
std::runtime_errorif R is not positive definite.

Based on the Schur Vector approach outlined in this paper: "On the Numerical Solution of the Discrete-Time Algebraic Riccati Equation" by Thrasyvoulos Pappas, Alan J. Laub, and Nils R. Sandell, in TAC, 1980, http://ieeexplore.ieee.org/stamp/stamp.jsp?arnumber=1102434

Note: When, for example, n = 100, m = 80, and entries of A, B, Q_half, R_half are sampled from standard normal distributions, where Q = Q_half'*Q_half and similar for R, the absolute error of the solution is 10^{-6}, while the absolute error of the solution computed by Matlab is 10^{-8}.

TODO(weiqiao.han): I may overwrite the RealQZ function to improve the accuracy, together with more thorough tests.

\[ A'XA - X - A'XB(B'XB+R)^{-1}B'XA + Q = 0 \]

Exceptions
std::runtime_errorif Q is not positive semi-definite.
std::runtime_errorif R is not positive definite.

Based on the Schur Vector approach outlined in this paper: "On the Numerical Solution of the Discrete-Time Algebraic Riccati Equation" by Thrasyvoulos Pappas, Alan J. Laub, and Nils R. Sandell

Here is the call graph for this function:

Here is the caller graph for this function:

drake::math::Gradient<Eigen::Matrix<typename Derived::Scalar, 3, 3>, drake::kQuaternionSize>::type drake::math::dquat2rotmat ( const Eigen::MatrixBase< Derived > &  quaternion)

Computes the gradient of the function that converts a unit length quaternion to a rotation matrix.

Parameters
quaternionA unit length quaternion [w;x;y;z]
Returns
The gradient

Here is the call graph for this function:

drake::math::DRAKE_DEPRECATED ( "This code is deprecated per issue #8323. ""Use RotationMatrix::ToQuaternion."  )

(Deprecated), use math::RotationMatrix::ToQuaternion

Here is the call graph for this function:

drake::math::DRAKE_DEPRECATED ( "This code is deprecated. Use Eigen's Quaternion constructor ""Quaternion(w, x, y, z) -- not a home-brew 4-element vectors"  )

Here is the caller graph for this function:

drake::math::DRAKE_DEPRECATED ( "This code is deprecated per issue #8323. ""Use math::RotationMatrix(R) or math::rotmat2rpy(R) ""or math::RotationMatrix::ToQuaternionAsVector4."  )
drake::math::DRAKE_DEPRECATED ( "This code is deprecated per issue #8323. ""Use constructor RollPitchYaw(quaternion)."  ) const

(Deprecated), use math::RollPitchYaw(quaternion).

Here is the call graph for this function:

DRAKE_DEPRECATED ( "This code is deprecated per issue #8323. ""Use RollPitchYaw::ToQuaternion()."  )

(Deprecated), use math::RollPitchYaw(rpy).ToQuaternion().

Here is the call graph for this function:

drake::math::Gradient< Eigen::Matrix<typename DerivedR::Scalar, drake::kQuaternionSize, 1>, DerivedDR::ColsAtCompileTime>::type drake::math::drotmat2quat ( const Eigen::MatrixBase< DerivedR > &  R,
const Eigen::MatrixBase< DerivedDR > &  dR 
)

Computes the gradient of the function that converts rotation matrix to quaternion.

Parameters
RA 3 x 3 rotation matrix
dRA 9 x N matrix, dR(i,j) is the gradient of R(i) w.r.t x_var(j)
Returns
The gradient G. G is a 4 x N matrix G(0,j) is the gradient of w w.r.t x_var(j) G(1,j) is the gradient of x w.r.t x_var(j) G(2,j) is the gradient of y w.r.t x_var(j) G(3,j) is the gradient of z w.r.t x_var(j)

Here is the call graph for this function:

drake::math::Gradient< Eigen::Matrix<typename DerivedR::Scalar, drake::kRpySize, 1>, DerivedDR::ColsAtCompileTime>::type drake::math::drotmat2rpy ( const Eigen::MatrixBase< DerivedR > &  R,
const Eigen::MatrixBase< DerivedDR > &  dR 
)

Computes the gradient of the function that converts a rotation matrix to body-fixed z-y'-x'' Euler angles.

Parameters
RA 3 x 3 rotation matrix
dRA 9 x N matrix, dR(i,j) is the gradient of R(i) w.r.t x(j)
Returns
The gradient G. G is a 3 x N matrix. G(0,j) is the gradient of roll w.r.t x(j) G(1,j) is the gradient of pitch w.r.t x(j) G(2,j) is the gradient of yaw w.r.t x(j)

Here is the call graph for this function:

Eigen::Matrix<typename Derived::Scalar, 4, 1> drake::math::expmap2quat ( const Eigen::MatrixBase< Derived > &  v)

Here is the call graph for this function:

Here is the caller graph for this function:

Eigen::Matrix<typename Derived::Scalar, Eigen::Dynamic, Eigen::Dynamic> drake::math::getSubMatrixGradient ( const Eigen::MatrixBase< Derived > &  dM,
const std::vector< int > &  rows,
const std::vector< int > &  cols,
typename Derived::Index  M_rows,
int  q_start = 0,
typename Derived::Index  q_subvector_size = -1 
)
GetSubMatrixGradientArray<QSubvectorSize, Derived, NRows, NCols>::type drake::math::getSubMatrixGradient ( const Eigen::MatrixBase< Derived > &  dM,
const std::array< int, NRows > &  rows,
const std::array< int, NCols > &  cols,
typename Derived::Index  M_rows,
int  q_start = 0,
typename Derived::Index  q_subvector_size = QSubvectorSize 
)
GetSubMatrixGradientSingleElement<QSubvectorSize, Derived>::type drake::math::getSubMatrixGradient ( const Eigen::MatrixBase< Derived > &  dM,
int  row,
int  col,
typename Derived::Index  M_rows,
typename Derived::Index  q_start = 0,
typename Derived::Index  q_subvector_size = QSubvectorSize 
)
void drake::math::gradientMatrixToAutoDiff ( const Eigen::MatrixBase< DerivedGradient > &  gradient,
Eigen::MatrixBase< DerivedAutoDiff > &  auto_diff_matrix 
)

Here is the caller graph for this function:

int GrayCodeToInteger ( const Eigen::Ref< const Eigen::VectorXi > &  gray_code)

Converts the Gray code to an integer.

For example (0, 0) -> 0 (0, 1) -> 1 (1, 1) -> 2 (1, 0) -> 3

Parameters
gray_codeThe N-digit Gray code, where N is gray_code.rows()
Returns
The integer represented by the Gray code gray_code.

Here is the caller graph for this function:

decltype(auto) drake::math::hessian ( F &&  f,
Arg &&  x 
)

Computes a matrix of AutoDiffScalars from which the value, Jacobian, and Hessian of a function

\[ f:\mathbb{R}^{n\times m}\rightarrow\mathbb{R}^{p\times q} \]

(f: R^n*m -> R^p*q) can be extracted.

The output is a matrix of nested AutoDiffScalars, being the result of calling jacobian on a function that returns the output of jacobian, called on f.

MaxChunkSizeOuter and MaxChunkSizeInner can be used to control chunk sizes (see jacobian).

See jacobian for requirements on the function f and the argument x.

Parameters
ffunction
xfunction argument value at which Hessian will be evaluated
Returns
AutoDiffScalar matrix corresponding to the Hessian of f evaluated at x
void drake::math::initializeAutoDiff ( const Eigen::MatrixBase< Derived > &  val,
Eigen::MatrixBase< DerivedAutoDiff > &  auto_diff_matrix,
Eigen::DenseIndex  num_derivatives = Eigen::Dynamic,
Eigen::DenseIndex  deriv_num_start = 0 
)

Initialize a single autodiff matrix given the corresponding value matrix.

Set the values of auto_diff_matrix to be equal to val, and for each element i of auto_diff_matrix, resize the derivatives vector to num_derivatives, and set derivative number deriv_num_start + i to one (all other elements of the derivative vector set to zero).

Parameters
[in]mat'regular' matrix of values
[out]retAutoDiff matrix
[in]num_derivativesthe size of the derivatives vector
Default: the size of mat
[in]deriv_num_startstarting index into derivative vector (i.e. element deriv_num_start in derivative vector corresponds to mat(0, 0)).
Default: 0

Here is the caller graph for this function:

AutoDiffMatrixType<Derived, Nq> drake::math::initializeAutoDiff ( const Eigen::MatrixBase< Derived > &  mat,
Eigen::DenseIndex  num_derivatives = -1,
Eigen::DenseIndex  deriv_num_start = 0 
)

Initialize a single autodiff matrix given the corresponding value matrix.

Create autodiff matrix that matches mat in size with derivatives of compile time size Nq and runtime size num_derivatives. Set its values to be equal to val, and for each element i of auto_diff_matrix, set derivative number deriv_num_start + i to one (all other derivatives set to zero).

Parameters
[in]mat'regular' matrix of values
[in]num_derivativesthe size of the derivatives vector
Default: the size of mat
[in]deriv_num_startstarting index into derivative vector (i.e. element deriv_num_start in derivative vector corresponds to mat(0, 0)).
Default: 0
Returns
AutoDiff matrix

Here is the call graph for this function:

void drake::math::initializeAutoDiffGivenGradientMatrix ( const Eigen::MatrixBase< Derived > &  val,
const Eigen::MatrixBase< DerivedGradient > &  gradient,
Eigen::MatrixBase< DerivedAutoDiff > &  auto_diff_matrix 
)

Initializes an autodiff matrix given a matrix of values and gradient matrix.

Parameters
[in]valvalue matrix
[in]gradientgradient matrix; the derivatives of val(j) are stored in row j of the gradient matrix.
[out]autodiff_matrixmatrix of AutoDiffScalars with the same size as val

Here is the caller graph for this function:

AutoDiffMatrixType<Derived, DerivedGradient::ColsAtCompileTime> drake::math::initializeAutoDiffGivenGradientMatrix ( const Eigen::MatrixBase< Derived > &  val,
const Eigen::MatrixBase< DerivedGradient > &  gradient 
)

Creates and initializes an autodiff matrix given a matrix of values and gradient matrix.

Parameters
[in]valvalue matrix
[in]gradientgradient matrix; the derivatives of val(j) are stored in row j of the gradient matrix.
Returns
autodiff_matrix matrix of AutoDiffScalars with the same size as val

Here is the call graph for this function:

std::tuple<AutoDiffMatrixType< Deriveds, internal::totalSizeAtCompileTime<Deriveds...>)>...> drake::math::initializeAutoDiffTuple ( const Eigen::MatrixBase< Deriveds > &...  args)

Given a series of Eigen matrices, create a tuple of corresponding AutoDiff matrices with values equal to the input matrices and properly initialized derivative vectors.

The size of the derivative vector of each element of the matrices in the output tuple will be the same, and will equal the sum of the number of elements of the matrices in args. If all of the matrices in args have fixed size, then the derivative vectors will also have fixed size (being the sum of the sizes at compile time of all of the input arguments), otherwise the derivative vectors will have dynamic size. The 0th element of the derivative vectors will correspond to the derivative with respect to the 0th element of the first argument. Subsequent derivative vector elements correspond first to subsequent elements of the first input argument (traversed first by row, then by column), and so on for subsequent arguments.

Parameters
argsa series of Eigen matrices
Returns
a tuple of properly initialized AutoDiff matrices corresponding to args

Here is the call graph for this function:

Here is the caller graph for this function:

std::array<int, Size> drake::math::intRange ( int  start)
bool drake::math::is_quaternion_in_canonical_form ( const Eigen::Quaternion< T > &  quat)

This function tests whether a quaternion is in "canonical form" meaning that it tests whether the quaternion [w, x, y, z] has a non-negative w value.

Example: [-0.3, +0.4, +0.5, +0.707] is not in canonical form. Example: [+0.3, -0.4, -0.5, -0.707] is in canonical form.

Parameters
quatQuaternion [w, x, y, z] that relates two right-handed orthogonal unitary bases e.g., Ax, Ay, Az (A) to Bx, By, Bz (B). Note: quat is analogous to the rotation matrix R_AB.
Returns
true if quat.w() is nonnegative (in canonical form), else false.

Here is the caller graph for this function:

bool drake::math::IsBothQuaternionAndQuaternionDtOK ( const Eigen::Quaternion< T > &  quat,
const Vector4< T > &  quatDt,
const double  tolerance 
)

This function tests if a quaternion satisfies the time-derivative constraint specified in [Kane, 1983] Section 1.13, equation 13, page 59.

A quaternion [w, x, y, z] must satisfy w^2 + x^2 + y^2 + z^2 = 1, hence its time-derivative must satisfy 2*(w*ẇ + x*ẋ + y*ẏ + z*ż) = 0. Note: To accurately test whether the time-derivative quaternion constraint is satisfied, the quaternion constraint is also tested to be accurate.

Parameters
quatQuaternion [w, x, y, z] that relates two right-handed orthogonal unitary bases e.g., Ax, Ay, Az (A) to Bx, By, Bz (B). Note: A quaternion like quat_AB is analogous to the rotation matrix R_AB.
quatDtTime-derivative of quat, i.e., [ẇ, ẋ, ẏ, ż].
toleranceTolerance for quaternion constraints.
Returns
true if both of the two previous constraints are within tolerance.

Here is the call graph for this function:

Here is the caller graph for this function:

bool drake::math::IsQuaternionAndQuaternionDtEqualAngularVelocityExpressedInB ( const Eigen::Quaternion< T > &  quat,
const Vector4< T > &  quatDt,
const Vector3< T > &  w_B,
const double  tolerance 
)

This function tests if a quaternion and a quaternions time-derivative can calculate and match an angular velocity to within a tolerance.

Note: This function first tests if the quaternion [w, x, y, z] satisifies w^2 + x^2 + y^2 + z^2 = 1 (to within tolerance) and if its time-derivative satisfies w*ẇ + x*ẋ + y*ẏ + z*ż = 0 (to within tolerance). Lastly, it tests if each element of the angular velocity calculated from quat and quatDt is within tolerance of w_B (described below).

Parameters
quatQuaternion [w, x, y, z] that relates two right-handed orthogonal unitary bases e.g., Ax, Ay, Az (A) to Bx, By, Bz (B). Note: A quaternion like quat_AB is analogous to the rotation matrix R_AB.
quatDtTime-derivative of quat, i.e., [ẇ, ẋ, ẏ, ż].
w_BRigid body B's angular velocity in frame A, expressed in B.
toleranceTolerance for quaternion constraints.
Returns
true if all three of the previous constraints are within tolerance.

Here is the call graph for this function:

Here is the caller graph for this function:

bool drake::math::IsQuaternionValid ( const Eigen::Quaternion< T > &  quat,
const double  tolerance 
)

This function tests if a quaternion satisfies the quaternion constraint specified in [Kane, 1983] Section 1.3, equation 4, page 12, i.e., a quaternion [w, x, y, z] must satisfy: w^2 + x^2 + y^2 + z^2 = 1.

Parameters
quatQuaternion [w, x, y, z] that relates two right-handed orthogonal unitary bases e.g., Ax, Ay, Az (A) to Bx, By, Bz (B). Note: A quaternion like quat_AB is analogous to the rotation matrix R_AB.
toleranceTolerance for quaternion constraint, i.e., how much is w^2 + x^2 + y^2 + z^2 allowed to differ from 1.
Returns
true if the quaternion constraint is satisfied within tolerance.

Here is the call graph for this function:

Here is the caller graph for this function:

bool drake::math::IsSymmetric ( const Eigen::MatrixBase< Derived > &  matrix)

Determines if a matrix is symmetric.

If std::equal_to<>()(matrix(i, j), matrix(j, i)) is true for all i, j, then the matrix is symmetric.

Here is the caller graph for this function:

bool drake::math::IsSymmetric ( const Eigen::MatrixBase< Derived > &  matrix,
const typename Derived::Scalar &  precision 
)

Determines if a matrix is symmetric based on whether the difference between matrix(i, j) and matrix(j, i) is smaller than precision for all i, j.

The precision is absolute. Matrix with nan or inf entries is not allowed.

decltype(auto) drake::math::jacobian ( F &&  f,
Arg &&  x 
)

Computes a matrix of AutoDiffScalars from which both the value and the Jacobian of a function

\[ f:\mathbb{R}^{n\times m}\rightarrow\mathbb{R}^{p\times q} \]

(f: R^n*m -> R^p*q) can be extracted.

The derivative vector for each AutoDiffScalar in the output contains the derivatives with respect to all components of the argument \( x \).

The return type of this function is a matrix with the `best' possible AutoDiffScalar scalar type, in the following sense:

  • If the number of derivatives can be determined at compile time, the AutoDiffScalar derivative vector will have that fixed size.
  • If the maximum number of derivatives can be determined at compile time, the AutoDiffScalar derivative vector will have that maximum fixed size.
  • If neither the number, nor the maximum number of derivatives can be determined at compile time, the output AutoDiffScalar derivative vector will be dynamically sized.

f should have a templated call operator that maps an Eigen matrix argument to another Eigen matrix. The scalar type of the output of \( f \) need not match the scalar type of the input (useful in recursive calls to the function to determine higher order derivatives). The easiest way to create an f is using a C++14 generic lambda.

The algorithm computes the Jacobian in chunks of up to MaxChunkSize derivatives at a time. This has three purposes:

  • It makes it so that derivative vectors can be allocated on the stack, eliminating dynamic allocations and improving performance if the maximum number of derivatives cannot be determined at compile time.
  • It gives control over, and limits the number of required instantiations of the call operator of f and all the functions it calls.
  • Excessively large derivative vectors can result in CPU capacity cache misses; even if the number of derivatives is fixed at compile time, it may be better to break up into chunks if that means that capacity cache misses can be prevented.
Parameters
ffunction
xfunction argument value at which Jacobian will be evaluated
Returns
AutoDiffScalar matrix corresponding to the Jacobian of f evaluated at x.

Here is the call graph for this function:

MatGradMult<DerivedDA, DerivedB>::type drake::math::matGradMult ( const Eigen::MatrixBase< DerivedDA > &  dA,
const Eigen::MatrixBase< DerivedB > &  B 
)

Here is the caller graph for this function:

MatGradMultMat<DerivedA, DerivedB, DerivedDA>::type drake::math::matGradMultMat ( const Eigen::MatrixBase< DerivedA > &  A,
const Eigen::MatrixBase< DerivedB > &  B,
const Eigen::MatrixBase< DerivedDA > &  dA,
const Eigen::MatrixBase< DerivedDB > &  dB 
)

Here is the caller graph for this function:

void drake::math::NormalizeVector ( const Eigen::MatrixBase< Derived > &  x,
typename Derived::PlainObject &  x_norm,
typename drake::math::Gradient< Derived, Derived::RowsAtCompileTime, 1 >::type *  dx_norm = nullptr,
typename drake::math::Gradient< Derived, Derived::RowsAtCompileTime, 2 >::type *  ddx_norm = nullptr 
)

Computes the normalized vector, optionally with its gradient and second derivative.

Parameters
[in]xAn N x 1 vector to be normalized. Must not be zero.
[out]x_normThe normalized vector (N x 1).
[out]dx_normIf non-null, returned as an N x N matrix, where dx_norm(i,j) = D x_norm(i)/D x(j).
[out]ddx_normIf non-null, and dx_norm is non-null, returned as an N^2 x N matrix, where ddx_norm.col(j) = D dx_norm/D x(j), with dx_norm stacked columnwise.

(D x / D y above means partial derivative of x with respect to y.)

Here is the call graph for this function:

Here is the caller graph for this function:

Matrix3<typename Derived::Scalar> drake::math::ProjectMatToOrthonormalMat ( const Eigen::MatrixBase< Derived > &  M)

(Deprecated), use math::RotationMatrix::ProjectToRotationMatrix

Here is the call graph for this function:

Matrix3<typename Derived::Scalar> drake::math::ProjectMatToRotMat ( const Eigen::MatrixBase< Derived > &  M)

(Deprecated), use math::RotationMatrix::ProjectToRotationMatrix

Here is the call graph for this function:

double drake::math::ProjectMatToRotMatWithAxis ( const Eigen::MatrixBase< Derived > &  M,
const Eigen::Ref< const Eigen::Vector3d > &  axis,
double  angle_lb,
double  angle_ub 
)

Projects a 3 x 3 matrix M onto SO(3).

The projected rotation matrix R has a given rotation axis a, and its rotation angle θ is bounded as angle_lb <= θ <= angle_ub. One use case for this function is to reconstruct the rotation matrix for a revolute joint with joint limits.

See also
GlobalInverseKinematics for an usage of this function. We can formulate this as an optimization problem
  min_θ trace((R - M)ᵀ*(R - M))
  subject to R = I + sinθ * A + (1 - cosθ) * A²   (1)
             angle_lb <= θ <= angle_ub
where A is the cross product matrix of the rotation axis a.
  A = [ 0  -a₃  a₂]
      [ a₃  0  -a₁]
      [-a₂  a₁  0 ]
Equation (1) is the Rodriguez Formula, to compute the rotation matrix from the rotation axis a and the rotation angle θ. For more details, refer to http://mathworld.wolfram.com/RodriguesRotationFormula.html The objective function can be simplified as
  max_θ trace(Rᵀ * M + Mᵀ * R)
By substituting the matrix R with the axis-angle representation, the optimization problem is formulated as
  max_θ sinθ * trace(Aᵀ*M) - cosθ * trace(Mᵀ * A²)
  subject to angle_lb <= θ <= angle_ub
By introducing α = atan2(-trace(Mᵀ * A²), trace(Aᵀ*M)), we can compute the optimal θ as
θ = π/2 + 2kπ - α, if angle_lb <= π/2 + 2kπ - α <= angle_ub, k ∈ ℤ
else
θ = angle_lb if sin(angle_lb + α) >= sin(angle_ub + α)
θ = angle_ub if sin(angle_lb + α) < sin(angle_ub + α)
Template Parameters
DerivedA 3 x 3 matrix
Parameters
MThe matrix to be projected.
axisThe axis of the rotation matrix. A unit length vector.
angle_lbThe lower bound of the rotation angle.
angle_ubThe upper bound of the rotation angle.
Returns
The rotation angle of the projected matrix.
Precondition
angle_ub >= angle_lb. Throw std::runtime_error if these bounds are violated.

Here is the call graph for this function:

Here is the caller graph for this function:

Eigen::Matrix<typename DerivedQ::Scalar, 3, 1> drake::math::quat2expmap ( const Eigen::MatrixBase< DerivedQ > &  q)

Here is the call graph for this function:

Here is the caller graph for this function:

void drake::math::quat2expmapSequence ( const Eigen::MatrixBase< DerivedQ > &  quat,
const Eigen::MatrixBase< DerivedQ > &  quat_dot,
Eigen::MatrixBase< DerivedE > &  expmap,
Eigen::MatrixBase< DerivedE > &  expmap_dot 
)

Here is the call graph for this function:

Matrix3<typename Derived::Scalar> drake::math::quat2rotmat ( const Eigen::MatrixBase< Derived > &  v)

Computes the rotation matrix from quaternion representation.

Template Parameters
DerivedAn Eigen derived type, e.g., an Eigen Vector3d.
Parameters
quaternion4 x 1 unit length quaternion, q=[w;x;y;z]
Returns
3 x 3 rotation matrix (Deprecated), use math::RotationMatrix(quaternion).

Here is the call graph for this function:

Here is the caller graph for this function:

Vector4<typename Derived::Scalar> drake::math::quatConjugate ( const Eigen::MatrixBase< Derived > &  q)

Here is the caller graph for this function:

Vector4<typename Derived1::Scalar> drake::math::quatDiff ( const Eigen::MatrixBase< Derived1 > &  q1,
const Eigen::MatrixBase< Derived2 > &  q2 
)

Here is the call graph for this function:

Here is the caller graph for this function:

Derived1::Scalar drake::math::quatDiffAxisInvar ( const Eigen::MatrixBase< Derived1 > &  q1,
const Eigen::MatrixBase< Derived2 > &  q2,
const Eigen::MatrixBase< DerivedU > &  u 
)

Here is the call graph for this function:

Here is the caller graph for this function:

Eigen::Quaternion<T> drake::math::QuaternionToCanonicalForm ( const Eigen::Quaternion< T > &  quat)

This function returns a quaternion in its "canonical form" meaning that it returns a quaternion [w, x, y, z] with a non-negative w.

For example, if passed a quaternion [-0.3, +0.4, +0.5, +0.707], the function returns the quaternion's canonical form [+0.3, -0.4, -0.5, -0.707].

Parameters
quatQuaternion [w, x, y, z] that relates two right-handed orthogonal unitary bases e.g., Ax, Ay, Az (A) to Bx, By, Bz (B). Note: quat is analogous to the rotation matrix R_AB.
Returns
Canonical form of quat, which means that either the original quat is returned or a quaternion representing the same orientation but with negated [w, x, y, z], to ensure a positive w in returned quaternion.

Here is the call graph for this function:

Here is the caller graph for this function:

Derived::Scalar drake::math::quatNorm ( const Eigen::MatrixBase< Derived > &  q)

Here is the call graph for this function:

Vector4<typename Derived1::Scalar> drake::math::quatProduct ( const Eigen::MatrixBase< Derived1 > &  q1,
const Eigen::MatrixBase< Derived2 > &  q2 
)

Here is the caller graph for this function:

Vector3<typename DerivedV::Scalar> drake::math::quatRotateVec ( const Eigen::MatrixBase< DerivedQ > &  q,
const Eigen::MatrixBase< DerivedV > &  v 
)

Here is the call graph for this function:

Here is the caller graph for this function:

void drake::math::resizeDerivativesToMatchScalar ( Eigen::MatrixBase< Derived > &  mat,
const typename Derived::Scalar &  scalar 
)

Resize derivatives vector of each element of a matrix to to match the size of the derivatives vector of a given scalar.

If the mat and scalar inputs are AutoDiffScalars, resize the derivatives vector of each element of the matrix mat to match the number of derivatives of the scalar. This is useful in functions that return matrices that do not depend on an AutoDiffScalar argument (e.g. a function with a constant output), while it is desired that information about the number of derivatives is preserved.

Parameters
matmatrix, for which the derivative vectors of the elements will be resized
scalarscalar to match the derivative size vector against.

Here is the caller graph for this function:

Vector3<typename Derived::Scalar> drake::math::rotmat2rpy ( const Eigen::MatrixBase< Derived > &  R)

Computes SpaceXYZ Euler angles from rotation matrix.

Template Parameters
DerivedAn Eigen derived type, e.g., an Eigen Vector3d.
Parameters
R3x3 rotation matrix.
Returns
3x1 SpaceXYZ Euler angles (called roll-pitch-yaw by ROS). Note: SpaceXYZ roll-pitch-yaw is equivalent to BodyZYX yaw-pitch-roll. http://answers.ros.org/question/58863/incorrect-rollpitch-yaw-values-using-getrpy/
See also
rpy2rotmat

Here is the call graph for this function:

Here is the caller graph for this function:

Matrix3<typename Derived::Scalar> drake::math::rpy2rotmat ( const Eigen::MatrixBase< Derived > &  rpy)

Here is the caller graph for this function:

T1 drake::math::saturate ( const T1 &  value,
const T2 &  low,
const T3 &  high 
)

Saturates the input value between upper and lower bounds.

If value is within [low, high] then return it; else return the boundary.

Here is the call graph for this function:

Here is the caller graph for this function:

void drake::math::setSubMatrixGradient ( Eigen::MatrixBase< DerivedA > &  dM,
const Eigen::MatrixBase< DerivedB > &  dM_submatrix,
const std::vector< int > &  rows,
const std::vector< int > &  cols,
typename DerivedA::Index  M_rows,
typename DerivedA::Index  q_start = 0,
typename DerivedA::Index  q_subvector_size = -1 
)
void drake::math::setSubMatrixGradient ( Eigen::MatrixBase< DerivedA > &  dM,
const Eigen::MatrixBase< DerivedB > &  dM_submatrix,
const std::array< int, NRows > &  rows,
const std::array< int, NCols > &  cols,
typename DerivedA::Index  M_rows,
typename DerivedA::Index  q_start = 0,
typename DerivedA::Index  q_subvector_size = QSubvectorSize 
)
void drake::math::setSubMatrixGradient ( Eigen::MatrixBase< DerivedDM > &  dM,
const Eigen::MatrixBase< DerivedDMSub > &  dM_submatrix,
int  row,
int  col,
typename DerivedDM::Index  M_rows,
typename DerivedDM::Index  q_start = 0,
typename DerivedDM::Index  q_subvector_size = QSubvectorSize 
)
Vector4<Scalar> drake::math::Slerp ( const Eigen::MatrixBase< Derived1 > &  q1,
const Eigen::MatrixBase< Derived2 > &  q2,
const Scalar &  interpolation_parameter 
)

Q = Slerp(q1, q2, f) Spherical linear interpolation between two quaternions This function uses the implementation given in Algorithm 8 of [1].

Parameters
q1Initial quaternion (w, x, y, z)
q2Final quaternion (w, x, y, z)
interpolation_parameterbetween 0 and 1 (inclusive)
Return values
QInterpolated quaternion(s). 4-by-1 vector.

[1] Kuffner, J.J., "Effective sampling and distance metrics for 3D rigid body path planning," Robotics and Automation, 2004. Proceedings. ICRA '04. 2004 IEEE International Conference on , vol.4, no., pp.3993, 3998 Vol.4, April 26-May 1, 2004 doi: 10.1109/ROBOT.2004.1308895

Here is the call graph for this function:

Here is the caller graph for this function:

void drake::math::SparseMatrixToRowColumnValueVectors ( const Derived &  matrix,
std::vector< Eigen::Index > &  row_indices,
std::vector< Eigen::Index > &  col_indices,
std::vector< typename Derived::Scalar > &  val 
)

For a sparse matrix, return the row indices, the column indices, and value of the non-zero entries.

For example, the matrix

\[ mat = \begin{bmatrix} 1 & 0 & 2\\ 0 & 3 & 4\end{bmatrix} \]

has

\[ row = \begin{bmatrix} 0 & 1 & 0 & 1\end{bmatrix}\\ col = \begin{bmatrix} 0 & 1 & 2 & 2\end{bmatrix}\\ val = \begin{bmatrix} 1 & 3 & 2 & 4\end{bmatrix} \]

Parameters
[in]matrixthe input sparse matrix
[out]row_indicesa vector containing the row indices of the non-zero entries
[out]col_indicesa vector containing the column indices of the non-zero entries
[out]vala vector containing the values of the non-zero entries.
std::vector<Eigen::Triplet<typename Derived::Scalar> > drake::math::SparseMatrixToTriplets ( const Derived &  matrix)

For a sparse matrix, return a vector of triplets, such that we can reconstruct the matrix using setFromTriplet function.

Parameters
matrixA sparse matrix
Returns
A triplet with the row, column and value of the non-zero entries. See https://eigen.tuxfamily.org/dox/group__TutorialSparse.html for more information on the triplet
drake::MatrixX<typename Derived::Scalar> drake::math::ToSymmetricMatrixFromLowerTriangularColumns ( const Eigen::MatrixBase< Derived > &  lower_triangular_columns)

Given a column vector containing the stacked columns of the lower triangular part of a square matrix, returning a symmetric matrix whose lower triangular part is the same as the original matrix.

Here is the call graph for this function:

Here is the caller graph for this function:

Eigen::Matrix<typename Derived::Scalar, rows, rows> drake::math::ToSymmetricMatrixFromLowerTriangularColumns ( const Eigen::MatrixBase< Derived > &  lower_triangular_columns)

Given a column vector containing the stacked columns of the lower triangular part of a square matrix, returning a symmetric matrix whose lower triangular part is the same as the original matrix.

Template Parameters
rowsThe number of rows in the symmetric matrix.

Here is the call graph for this function:

Derived::PlainObject drake::math::transposeGrad ( const Eigen::MatrixBase< Derived > &  dX,
typename Derived::Index  rows_X 
)

Here is the caller graph for this function:

Eigen::AngleAxisd drake::math::UniformlyRandomAngleAxis ( Generator *  generator)

Generates a rotation (in the axis-angle representation) that rotates a point on the unit sphere to another point on the unit sphere with a uniform distribution over the sphere.

Justification for the algorithm can be found in, e.g.: Mervin E. Muller. 1959. A note on a method for generating points uniformly on n-dimensional spheres. Commun. ACM 2, 4 (April 1959), 19-20. DOI=http://dx.doi.org/10.1145/377939.377946

Here is the caller graph for this function:

Eigen::Quaterniond drake::math::UniformlyRandomQuaternion ( Generator *  generator)

Generates a rotation (in the quaternion representation) that rotates a point on the unit sphere to another point on the unit sphere with a uniform distribution over the sphere.

Here is the call graph for this function:

Here is the caller graph for this function:

RotationMatrix<double> drake::math::UniformlyRandomRotationMatrix ( Generator *  generator)

Generates a rotation (in the rotation matrix representation) that rotates a point on the unit sphere to another point on the unit sphere with a uniform distribution over the sphere.

Here is the call graph for this function:

Eigen::Vector3d drake::math::UniformlyRandomRPY ( Generator *  generator)

Generates a rotation (in the roll-pitch-yaw representation) that rotates a point on the unit sphere to another point on the unit sphere with a uniform distribution over the sphere.

Here is the call graph for this function:

Here is the caller graph for this function:

Eigen::Matrix3Xd UniformPtsOnSphereFibonacci ( int  num_points)

Deterministically generates approximate evenly distributed points on a unit sphere.

This method uses Fibonacci number. For the detailed math, please refer to http://stackoverflow.com/questions/9600801/evenly-distributing-n-points-on-a-sphere This algorithm generates the points in O(n) time, where n is the number of points.

Parameters
num_pointsThe number of points we want on the unit sphere.
Returns
The generated points.
Precondition
num_samples >= 1. Throw std::runtime_error if num_points < 1

Here is the call graph for this function:

drake::Matrix3<typename Derived::Scalar> drake::math::VectorToSkewSymmetric ( const Eigen::MatrixBase< Derived > &  p)

Here is the caller graph for this function:

T1 drake::math::wrap_to ( const T1 &  value,
const T2 &  low,
const T2 &  high 
)

For variables that are meant to be periodic, (e.g.

over a 2π interval), wraps value into the interval [low, high). Precisely, wrap_to returns: value + k*(high-low) for the unique integer value k that lands the output in the desired interval. low and high must be finite, and low < high.

Here is the call graph for this function:

Here is the caller graph for this function:

Matrix3<T> drake::math::XRotation ( const T &  theta)

(Deprecated), use math::RotationMatrix::MakeXRotation

Here is the call graph for this function:

Matrix3<T> drake::math::YRotation ( const T &  theta)

(Deprecated), use math::RotationMatrix::MakeYRotation

Here is the call graph for this function:

Matrix3<T> drake::math::ZRotation ( const T &  theta)

(Deprecated), use math::RotationMatrix::MakeZRotation

Here is the call graph for this function: