htransform.h File Reference

Detailed Description

Module to manage homogeneous transforms in R^3 used to manipulated 3D geometry.

See Also
htransform.c, THTransform, Tplot3d, plot3d.c.

Definition in file htransform.h.

Macros

#define DIM_SP   3
 Dimensions of R^3. More...
 
#define TX   0
 One of the types of homogeneous transforms in R^3. More...
 
#define TY   1
 One of the types of homogeneous transforms in R^3. More...
 
#define TZ   2
 One of the types of homogeneous transforms in R^3. More...
 
#define RX   3
 One of the types of homogeneous transforms in R^3. More...
 
#define RY   4
 One of the types of homogeneous transforms in R^3. More...
 
#define RZ   5
 One of the types of homogeneous transforms in R^3. More...
 
#define AXIS_X   0
 One of the dimension of R^3. More...
 
#define AXIS_Y   1
 One of the dimension of R^3. More...
 
#define AXIS_Z   2
 One of the dimension of R^3. More...
 
#define AXIS_H   3
 The homogeneous dimension in R^3. More...
 

Typedefs

typedef double THTransform [DIM_SP+1][DIM_SP+1]
 A homgeneous transform in R^3. More...
 

Functions

void HTransformIdentity (THTransform *t)
 Constructor. More...
 
void HTransformZero (THTransform *t)
 Constructor. More...
 
void HTransformCopy (THTransform *t_dst, THTransform *t_org)
 Copy constructor. More...
 
boolean HTransformIsIdentity (THTransform *t)
 Identify the identity matrix. More...
 
boolean HTransformIsTranslation (THTransform *t)
 Identify the translation matrices. More...
 
void HTransformTx (double tx, THTransform *t)
 Constructor. More...
 
void HTransformTy (double ty, THTransform *t)
 Constructor. More...
 
void HTransformTz (double tz, THTransform *t)
 Constructor. More...
 
void HTransformTxyz (double tx, double ty, double tz, THTransform *t)
 Constructor. More...
 
void HTransformRx (double rx, THTransform *t)
 Constructor. More...
 
void HTransformRy (double ry, THTransform *t)
 Constructor. More...
 
void HTransformRz (double rz, THTransform *t)
 Constructor. More...
 
void HTransformRx2 (double s, double c, THTransform *t)
 Constructor. More...
 
void HTransformRy2 (double s, double c, THTransform *t)
 Constructor. More...
 
void HTransformRz2 (double s, double c, THTransform *t)
 Constructor. More...
 
void HTransformScale (double s, THTransform *t)
 Constructor. More...
 
void HTransformScaleX (double s, THTransform *t)
 Constructor. More...
 
void HTransformScaleY (double s, THTransform *t)
 Constructor. More...
 
void HTransformScaleZ (double s, THTransform *t)
 Constructor. More...
 
void HTransformCreate (unsigned int dof_r3, double v, THTransform *t)
 Constructor. More...
 
void HTransformSetElement (unsigned int i, unsigned int j, double v, THTransform *t)
 Sets an element in a homogeneous transform. More...
 
double HTransformGetElement (unsigned int i, unsigned int j, THTransform *t)
 Gets an element in a homogeneous transform. More...
 
void HTransformFromVectors (double *x, double *y, double *z, double *h, THTransform *t)
 Defines a homogeneous transform from 4 vectors. More...
 
void HTransform2GLMatrix (double *m, THTransform *t)
 Defines a GL column-major matrix from a homogeneous transform. More...
 
void HTransformFromGLMatrix (double *m, THTransform *t)
 Defines homogeneous transform from a GL matrix. More...
 
void HTransformProduct (THTransform *t1, THTransform *t2, THTransform *t3)
 Product of two homogeneous transforms. More...
 
void HTransformAdd (THTransform *t1, THTransform *t2, THTransform *t3)
 Addition of two homogeneous transforms. More...
 
void HTransformInverse (THTransform *t, THTransform *ti)
 Inverse of a homogeneous transform. More...
 
void HTransformOrthonormalize (THTransform *t, THTransform *ta)
 Orthonormalizes the rotation part of a homogenouos transform. More...
 
void HTransformX2Vect (double sy, double sz, double *p1, double *p2, THTransform *t)
 Transform a unitary vector along the X axis to a generic vector. More...
 
void HTransformYawPitchRoll (double a, double b, double c, THTransform *t)
 Defines a rotation matrix from the yaw, pitch, and roll parameters. More...
 
void GetYawPitchRoll (double *a, double *b, double *c, THTransform *t)
 Recovers the Euler angles from a rotation matrix. More...
 
void HTransformTranspose (THTransform *t, THTransform *tt)
 Transpose of a homogeneous transform. More...
 
void HTransformAcumTrans (double tx, double ty, double tz, THTransform *t)
 Computes the result of multiplying a homogeneous transform by a translation matrix with parameters tx, ty, tz. More...
 
void HTransformAcumTrans2 (THTransform *t_in, double tx, double ty, double tz, THTransform *t)
 Computes the result of multiplying a homogeneous transform by a translation matrix with parameters tx, ty, tz. More...
 
void HTransformAcumRot (unsigned int type, double s, double c, THTransform *t)
 Computes the result of multiplying a homogeneous transform by a rotation matrix. More...
 
void HTransformAcumRot2 (THTransform *t_in, unsigned int type, double s, double c, THTransform *t)
 Computes the result of multiplying a homogeneous transform by a rotation matrix. More...
 
void HTransformApply (double *p_in, double *p_out, THTransform *t)
 Multiply a homogeneous transform and a vector. More...
 
void HTransformApplyRot (double *p_in, double *p_out, THTransform *t)
 Multiply the rotation part of the homogeneous transform and a vector. More...
 
void HTransformPrint (FILE *f, THTransform *t)
 Prints the a homogeneous transform to a file. More...
 
void HTransformPrintT (FILE *f, THTransform *t)
 Prints the transpose of a homogeneous transform to a file. More...
 
void HTransformPrettyPrint (FILE *f, THTransform *t)
 Prints a homogenoeus transform compactly. More...
 
void HTransformDelete (THTransform *t)
 Destructor. More...
 

Macro Definition Documentation

#define TX   0

One of the types of homogeneous transforms in R^3: Translation along X.

See Also
THTransform

Definition at line 35 of file htransform.h.

Referenced by DeriveTransSeq(), EvaluateTransSeq(), GetJointTransSeq(), HTransformCreate(), PrintTransSeq(), RecomputeScalarEquations(), SHTransformCreate(), SHTransformVarCreate(), and UpdateUsedDOF().

#define TY   1

One of the types of homogeneous transforms in R^3: Translation along Y.

See Also
THTransform

Definition at line 43 of file htransform.h.

Referenced by DeriveTransSeq(), EvaluateTransSeq(), HTransformCreate(), PrintTransSeq(), RecomputeScalarEquations(), SHTransformCreate(), SHTransformVarCreate(), and UpdateUsedDOF().

#define TZ   2

One of the types of homogeneous transforms in R^3: Translation along Z.

See Also
Thtransform

Definition at line 51 of file htransform.h.

Referenced by DeriveTransSeq(), EvaluateTransSeq(), HTransformCreate(), PrintTransSeq(), RecomputeScalarEquations(), SHTransformCreate(), SHTransformVarCreate(), and UpdateUsedDOF().

Typedef Documentation

typedef double THTransform[DIM_SP+1][DIM_SP+1]

A 4x4 matrix typically representing a homogenous matrix in R^3. This is also used to represent generic 4x4 matrices sometimes. Therefor operations do not assume any structure in the matrix.

See Also
htransform.h, htransform.c

Definition at line 129 of file htransform.h.

Function Documentation

void HTransformIdentity ( THTransform t)
void HTransformZero ( THTransform t)

Initializes a homogeneous transform with the zero matrix.

Parameters
tThe matrix to initialize.

Definition at line 75 of file htransform.c.

References MATRIX_INT_COPY.

Referenced by DeriveMEquation(), and DeriveTransSeq().

void HTransformCopy ( THTransform t_dst,
THTransform t_org 
)

Initializes a homogeneous transform from another transform.

Parameters
t_dstThe matrix to initialize.
t_orgThe matrix from where to copy.

Definition at line 83 of file htransform.c.

References MATRIX_INT_COPY.

Referenced by CDCallBackInfo(), CheckCollisionSolid(), CopyJoint(), CopyLink(), CopyMEquation(), CopyTrans(), GetJointDOFValues(), GetJointTransform(), InitCtTrans(), and NewFixJoint().

boolean HTransformIsIdentity ( THTransform t)

Identifies identity matrix (up to a tiny error).

Parameters
tThe matrix to check.
Returns
TRUE if the given matrix is the identity.

Definition at line 91 of file htransform.c.

References Distance(), and ZERO.

Referenced by AddCtTrans2TransSeq(), AddDispTrans2TransSeq(), AddPatchTrans2TransSeq(), AddTrans2TransSeq(), AddVarTrans2TransSeq(), HTransformPrettyPrint(), PrintMEquation(), PrintTransSeq(), SolidCorrection(), and TransformPolyhedron().

boolean HTransformIsTranslation ( THTransform t)

Identifies translation matrices (up to a tiny error), i.e., it does not include any rotation. It returns TRUE also for the identity (no rotation and no translation).

Parameters
tThe matrix to check.
Returns
TRUE if the given matrix includes no rotation.

Definition at line 96 of file htransform.c.

References AXIS_X, AXIS_Y, AXIS_Z, and ZERO.

Referenced by HasCtRotTransSeq(), and HTransformPrettyPrint().

void HTransformTx ( double  tx,
THTransform t 
)

Initializes a homogeneous transform as a translation along X.

Parameters
txThe displacement along X.
tThe matrix to initialize.

Definition at line 106 of file htransform.c.

References AXIS_H, AXIS_X, and MATRIX_INT_COPY.

Referenced by GetJointDOFValues(), GetJointTransform(), GetJointTransSeq(), HTransformCreate(), and main().

void HTransformTy ( double  ty,
THTransform t 
)

Initializes a homogeneous transform as a translation along Y.

Parameters
tyThe displacement along Y.
tThe matrix to initialize.

Definition at line 117 of file htransform.c.

References AXIS_H, AXIS_Y, and MATRIX_INT_COPY.

Referenced by HTransformCreate(), main(), and SolidCorrection().

void HTransformTz ( double  tz,
THTransform t 
)

Initializes a homogeneous transform as a translation along Z.

Parameters
tzThe displacement along Z.
tThe matrix to initialize.

Definition at line 128 of file htransform.c.

References AXIS_H, AXIS_Z, and MATRIX_INT_COPY.

Referenced by HTransformCreate(), and main().

void HTransformTxyz ( double  tx,
double  ty,
double  tz,
THTransform t 
)

Initializes a homogeneous transform as a translation.

Parameters
txThe displacement along X.
tyThe displacement along Y.
tzThe displacement along Z.
tThe matrix to initialize.

Definition at line 140 of file htransform.c.

References AXIS_H, AXIS_X, AXIS_Y, AXIS_Z, and MATRIX_INT_COPY.

Referenced by GenerateSphereOFF(), GetJointTransform(), HTransformX2Vect(), main(), NewPrismaticJoint(), NewSphericalJoint(), NewSphPrsSphJoint(), NewSphSphJoint(), and SolidCorrection().

void HTransformRx ( double  rx,
THTransform t 
)

Initializes a homogeneous transform as a rotation about X.

Parameters
rxThe rotation about X in radiants.
tThe matrix to initialize.

Definition at line 155 of file htransform.c.

References AXIS_Y, AXIS_Z, and MATRIX_INT_COPY.

Referenced by GetJointTransform(), HTransformCreate(), main(), NewRevoluteJoint(), and NewUniversalJoint().

void HTransformRy ( double  ry,
THTransform t 
)

Initializes a homogeneous transform as a rotation about Y.

Parameters
ryThe rotation about Y in radiants.
tThe matrix to initialize.

Definition at line 172 of file htransform.c.

References AXIS_X, AXIS_Z, and MATRIX_INT_COPY.

Referenced by GetJointTransform(), HTransformCreate(), HTransformX2Vect(), and main().

void HTransformRz ( double  rz,
THTransform t 
)

Initializes a homogeneous transform as a rotation about Z.

Parameters
rzThe rotation about Z in radiants.
tThe matrix to initialize.

Definition at line 189 of file htransform.c.

References AXIS_X, AXIS_Y, and MATRIX_INT_COPY.

Referenced by GenerateSphereOFF(), GetJointTransform(), HTransformCreate(), HTransformX2Vect(), HTransformYawPitchRoll(), and main().

void HTransformRx2 ( double  s,
double  c,
THTransform t 
)

Initializes a homogeneous transform as a rotation about X.

Parameters
sThe cosinus of the rotation about X.
cThe sinus of the rotation about X.
tThe matrix to initialize.

Definition at line 206 of file htransform.c.

References AXIS_Y, AXIS_Z, and MATRIX_INT_COPY.

Referenced by EvaluateTransSeq(), and GetJointTransSeq().

void HTransformRy2 ( double  s,
double  c,
THTransform t 
)

Initializes a homogeneous transform as a rotation about Y.

Parameters
sThe cosinus of the rotation about Y.
cThe sinus of the rotation about Y.
tThe matrix to initialize.

Definition at line 216 of file htransform.c.

References AXIS_X, AXIS_Z, and MATRIX_INT_COPY.

Referenced by EvaluateTransSeq().

void HTransformRz2 ( double  s,
double  c,
THTransform t 
)

Initializes a homogeneous transform as a rotation about Z.

Parameters
sThe cosinus of the rotation about Z.
cThe sinus of the rotation about Z.
tThe matrix to initialize.

Definition at line 226 of file htransform.c.

References AXIS_X, AXIS_Y, and MATRIX_INT_COPY.

Referenced by EvaluateTransSeq(), GetJointTransSeq(), and SolidCorrection().

void HTransformScale ( double  s,
THTransform t 
)

Initializes a homogeneous transform as a scale matrix.

Parameters
sThe ccale factor.
tThe matrix to initialize.

Definition at line 233 of file htransform.c.

References AXIS_X, AXIS_Y, AXIS_Z, and MATRIX_INT_COPY.

Referenced by main().

void HTransformScaleX ( double  s,
THTransform t 
)

Initializes a homogeneous transform as a scale matrix along the X axis.

Parameters
sThe ccale factor.
tThe matrix to initialize.

Definition at line 243 of file htransform.c.

References AXIS_X, and MATRIX_INT_COPY.

Referenced by main().

void HTransformScaleY ( double  s,
THTransform t 
)

Initializes a homogeneous transform as a scale matrix along the Y axis.

Parameters
sThe ccale factor.
tThe matrix to initialize.

Definition at line 251 of file htransform.c.

References AXIS_Y, and MATRIX_INT_COPY.

Referenced by main().

void HTransformScaleZ ( double  s,
THTransform t 
)

Initializes a homogeneous transform as a scale matrix along the Z axis.

Parameters
sThe ccale factor.
tThe matrix to initialize.

Definition at line 259 of file htransform.c.

References AXIS_Z, and MATRIX_INT_COPY.

Referenced by main().

void HTransformCreate ( unsigned int  dof_r3,
double  v,
THTransform t 
)

Initializes a homogeneous transform as translation/rotation in a given degree of freedom.

Parameters
dof_r3The degree of freedom: TX, TY, TZ, RX, RY, RZ.
vThe translation or rotation (in radiants) for the selected degree of freedom.
tThe matrix to initialize.

Definition at line 272 of file htransform.c.

References Error(), HTransformRx(), HTransformRy(), HTransformRz(), HTransformTx(), HTransformTy(), HTransformTz(), RX, RY, RZ, TX, TY, and TZ.

void HTransformSetElement ( unsigned int  i,
unsigned int  j,
double  v,
THTransform t 
)
inline

Sets an element in a homogeneous transform. Note that if the element value is not carefully computed the result of using this function could not be an homogeneous matrix anymore (recall that the 3x3 top left sub-matrix must be orthonormal).

Parameters
iRow of the element to set.
jColumn of the element to set.
vThe new value for the element
tThe matrix to update.

Definition at line 306 of file htransform.c.

References AXIS_H, and Error().

Referenced by ChangeLinkReferenceFrame(), DeriveTransSeq(), EvaluateTransSeq(), GetJointDOFValues(), GetTransform2LinkFLinks(), GetTransform2LinkLinks(), GetTransform2LinkQLinks(), HTransformX2Vect(), main(), and SHTransformEvaluate().

double HTransformGetElement ( unsigned int  i,
unsigned int  j,
THTransform t 
)
inline

Gets an element from a homogeneous transform.

Parameters
iRow of the element to get.
jColumn of the element to get.
tThe matrix to update.
Returns
The element at row i column j of homogeneous transform t.

Definition at line 323 of file htransform.c.

References AXIS_H, and Error().

Referenced by EvaluateMEquation(), GenerateJointEquations(), GenerateLinkSolutionFLinks(), GenerateLinkSolutionLinks(), GenerateLinkSolutionQLinks(), GetJointDOFValues(), NewFixJoint(), PQPTest(), PrintTransSeq(), SHTransformPostCtProduct(), SHTransformPreCtProduct(), and WorldDOF2Sol().

void HTransformFromVectors ( double *  x,
double *  y,
double *  z,
double *  h,
THTransform t 
)

Defines a homogeneous transform from 4 vectors defining its columns. No check is performed to verify if the rotation part (x,y,z columns) is actually a rotation matrix.

Parameters
xThe 3 values for the first column.
yThe 3 values for the second column.
zThe 3 values for the third column.
hThe 3 values for the fourth (translation) column.
tThe transform to define.

Definition at line 335 of file htransform.c.

References AXIS_H, AXIS_X, AXIS_Y, and AXIS_Z.

Referenced by Atoms2Transforms(), EvaluatePATrans(), EvaluateTransSeq(), GetJointDOFValues(), GetJointTransform(), NewRevoluteJoint(), and NewUniversalJoint().

void HTransform2GLMatrix ( double *  m,
THTransform t 
)

Defines a homogeneous transform a la GL (colum-major) from our format of homogeneous transforms.

Parameters
mSpace where to store the GL matrix.
tThe homogeneous transfrom.

Definition at line 351 of file htransform.c.

References AXIS_H.

Referenced by CheckCollisionSolid(), and InitSolidCD().

void HTransformFromGLMatrix ( double *  m,
THTransform t 
)

The reverse of HTransform2GLMatrix.

Parameters
mSpace from where to get the GL matrix.
tThe homogeneous transfrom to define.

Definition at line 373 of file htransform.c.

References AXIS_H.

void HTransformAdd ( THTransform t1,
THTransform t2,
THTransform t3 
)

Addition of two homogeneous transforms.

Parameters
t1The first transform to operate.
t2The second transform to operate.
t3The result.

Definition at line 437 of file htransform.c.

References DIM_SP, and MATRIX_INT_COPY.

void HTransformInverse ( THTransform t,
THTransform ti 
)
void HTransformOrthonormalize ( THTransform t,
THTransform ta 
)

Orthonormalizes the rotation part of a homogenouos transform.

This is to be used when we have an approximated homogeneous transform (for instance, a homogeneous transform that is derived after several floating point operations). If the rotation part is not orthonormal then the homogeneous transform is no longer a rigid transform (it produces deformations in the objects).

Right now this is implemented in a simple way

  • The first (column) vector of the matrix is normalized.
  • We substract the projectio of the first vector from the second one and we normalize it.
  • The third vector is the cross product of the two first vectors.

An alternative implementation could be based in a SVD of the rotation matrix: after the decomposition the diagonal matrix of (squared) eigenvalues is set to the identity and the matrix is reconstucted. In this way we could avoid the bias of adjusting the matrix from the first column.

Parameters
tThe transform to be orthonormalized
taThe resulting transform with the rotation part orthonormalized.

Definition at line 494 of file htransform.c.

References AXIS_H, and MATRIX_INT_COPY.

Referenced by GetTransform2LinkFLinks(), GetTransform2LinkLinks(), and GetTransform2LinkQLinks().

void HTransformX2Vect ( double  sy,
double  sz,
double *  p1,
double *  p2,
THTransform t 
)

Computes the homogeneous matrix that transform a unitary vector along the X axis into a generic vector.

The resulting transform includes, scaling, rotation and translation.

This function is typically used when plotting cylinders. This is the reason to include scale factors in Y and Z (to deform the unitary cylinder along X with radii 1 so that we get the cylinder in the desired position and with the desired witdh, typically the same in Y and Z). The length is controlled by the start-end points of the cilynder.

Parameters
syScale factor to apply in the Y axes.
szScale factor to apply in the Z axes.
p1Initial point of the vector.
p2Final point of the vector.
tThe resulting homogeneous transform.

Definition at line 539 of file htransform.c.

References Error(), HTransformIdentity(), HTransformProduct(), HTransformRy(), HTransformRz(), HTransformSetElement(), HTransformTxyz(), and M_PI_2.

Referenced by GenerateCylinderOFF(), GetJointDOFValues(), MoveJointFromTransforms(), NewInPatchJoint(), NewRevoluteJoint(), NewSphericalJoint(), NewUniversalJoint(), PlotCylinder(), and SolidCorrection().

void HTransformYawPitchRoll ( double  a,
double  b,
double  c,
THTransform t 
)

Defines de rotation matrix resulting from

R=Rz(a)Ry(b)Rx(c)

This is the invers of GetYawPitchRoll.

Parameters
aYaw. Rotation around z. This is the last being applied.
bPitch. Rotation around y. This is the second being applied.
cRoll. Rotation around x. This is the fisrt being applied.
tThe homogeneous transform with the rotation.

Definition at line 590 of file htransform.c.

References HTransformAcumRot(), HTransformRz(), RX, and RY.

Referenced by main().

void GetYawPitchRoll ( double *  a,
double *  b,
double *  c,
THTransform t 
)

Recovers the Yaw, Pitch, and Roll angles from the rotation part of a homogeneous matrix, i.e., recovers the parameters such that

   Rz(a)Ry(b)Rx(c)=R

with R the rotation part of the given transform. Note that in the expression

   Rz(a)Ry(b)Rx(c)

we first apply the roll, then the pitch, and finally the yaw.

The implementation uses the method taken from http://planning.cs.uiuc.edu/node103.html

This is the inverse of HTransformYawPitchRoll

Parameters
aOutput. Yaw. Rotation around z.
bOutput. Pitch. Rotation around y.
cOutput. Roll. Rotation around x.
tThe homogeneous transform with the rotation.

Definition at line 599 of file htransform.c.

References AXIS_X, AXIS_Y, and AXIS_Z.

Referenced by GetJointDOFValues(), and HTransformPrettyPrint().

void HTransformTranspose ( THTransform t,
THTransform tt 
)

Transpose of a homogeneous transform. Note that, in general the transpose of a homogeneous transform is not a homogeneous transform.

Parameters
tThe transform to transpose.
ttThe result with the transposed matrix.

Definition at line 618 of file htransform.c.

References DIM_SP, and MATRIX_INT_COPY.

Referenced by EvaluateTransSeq().

void HTransformAcumTrans ( double  tx,
double  ty,
double  tz,
THTransform t 
)

Computes the result of multiplying a homogeneous transform by a translation matrix with parameters tx, ty, tz. The result is accumulated in the given transform.

Parameters
txThe displacement along X.
tyThe displacement along Y.
tzThe displacement along Z.
tThe matrix where to accumulate the result of the product.
See Also
HTransformAcumTrans2

Definition at line 641 of file htransform.c.

References AXIS_H, AXIS_X, AXIS_Y, and AXIS_Z.

Referenced by EvaluateTransSeq(), HTransformAcumTrans2(), and MoveJointFromTransforms().

void HTransformAcumTrans2 ( THTransform t_in,
double  tx,
double  ty,
double  tz,
THTransform t 
)

Computes the result of multiplying a homogeneous transform by a translation matrix with parameters tx, ty, tz. The result is accumulated in the given transform.

Parameters
t_inThe base transform where to accumulate the translation.
txThe displacement along X.
tyThe displacement along Y.
tzThe displacement along Z.
tThe matrix where to store the result.
See Also
HTransformAcumTrans

Definition at line 656 of file htransform.c.

References HTransformAcumTrans(), and MATRIX_INT_COPY.

void HTransformAcumRot ( unsigned int  type,
double  s,
double  c,
THTransform t 
)

Computes the result of multiplying a homogeneous transform by a rotation matrix. The result is accumulated in the given transform.

Parameters
typeType of rotation to accumulate: RX, RY, RZ.
ssinus of the angle to rotate.
ccosinus of the angle to rotate.
tThe matrix from where to depart and where to store the result.
See Also
HTransformAcumRot2

Definition at line 673 of file htransform.c.

References AXIS_H, AXIS_X, AXIS_Y, AXIS_Z, Error(), RX, RY, and RZ.

Referenced by EvaluateTransSeq(), GetJointTransform(), HTransformAcumRot2(), and HTransformYawPitchRoll().

void HTransformAcumRot2 ( THTransform t_in,
unsigned int  type,
double  s,
double  c,
THTransform t 
)

Computes the result of multiplying a homogeneous transform by a rotation matrix.

Parameters
t_inThe base transform where to accumulate the translation.
typeType of rotation to accumulate: RX, RY, RZ.
ssinus of the angle to rotate.
ccosinus of the angle to rotate.
tThe matrix from where to store the result.
See Also
HTransformAcumRot

Definition at line 716 of file htransform.c.

References HTransformAcumRot(), and MATRIX_INT_COPY.

Referenced by MoveJointFromTransforms().

void HTransformApply ( double *  p_in,
double *  p_out,
THTransform t 
)
void HTransformApplyRot ( double *  p_in,
double *  p_out,
THTransform t 
)

Applies the rotation encoded in a homogeneous transform to a vector (i.e., rotates the vector).

Parameters
p_inA input vector in R^3.
p_outThe resulting vector in R^3.
tThe matrix with the rotation to apply to the input vector.

Definition at line 747 of file htransform.c.

References DIM_SP.

Referenced by EvaluateTransSeq(), GenerateJointSolution(), and PrintJointAxes().

void HTransformPrint ( FILE *  f,
THTransform t 
)

Prints the a homogeneous transform to a stream that can be stdtout.

Parameters
fThe file where to print.
tThe matrix to print.

Definition at line 768 of file htransform.c.

References DIM_SP.

Referenced by Atoms2Transforms(), and StoreCollisionInfoInt().

void HTransformPrintT ( FILE *  f,
THTransform t 
)

Prints the transpose of a homogeneous transform to a stream that can be stdtout. This is implemented because computer graphics community use transposed homogeneous transforms and, consequently, this is how geomview requires the transforms to be printed.

Parameters
fThe file where to print.
tThe matrix to print.
See Also
plot3d.h

Definition at line 794 of file htransform.c.

References DIM_SP.

Referenced by Move3dObject(), and PlotCylinder().

void HTransformPrettyPrint ( FILE *  f,
THTransform t 
)

Prints a homogenoeus transform compactly using basic transforms, whenever possible.

Parameters
fThe file where to print.
tThe matrix to print.

Definition at line 809 of file htransform.c.

References AXIS_H, AXIS_X, AXIS_Y, AXIS_Z, GetYawPitchRoll(), HTransformIsIdentity(), and HTransformIsTranslation().

Referenced by PrintJoint().