Drake
Drake C++ Documentation
Joint< T > Class Template Referenceabstract

Detailed Description

template<typename T>
class drake::multibody::Joint< T >

A Joint models the kinematical relationship which characterizes the possible relative motion between two bodies.

The two bodies connected by this Joint object are referred to as parent and child bodies. The parent/child ordering defines the sign conventions for the generalized coordinates and the coordinate ordering for multi-DOF joints. A Joint is a model of a physical kinematic constraint between two bodies, a constraint that in the real physical system does not specify a tree ordering.

BodyParentChildJoint.png

In Drake we define a frame F rigidly attached to the parent body P with pose X_PF and a frame M rigidly attached to the child body B with pose X_BM. A Joint object specifies a kinematic relation between frames F and M, which in turn imposes a kinematic relation between bodies P and B.

Typical joints include the ball joint, to allow unrestricted rotations about a given point, the revolute joint, that constraints two bodies to rotate about a given common axis, etc.

Consider the following example to build a simple pendulum system:

MultibodyPlant<double> plant(0.0);
// ... Code here to setup quantities below as mass, com, etc. ...
const RigidBody<double>& pendulum =
plant.AddRigidBody(SpatialInertia<double>(mass, com, unit_inertia));
// We will connect the pendulum body to the world using a RevoluteJoint.
// In this simple case the parent body P is the model's world body and frame
// F IS the world frame.
// Additionally, we need to specify the pose of frame M on the pendulum's
// body frame B.
// Say we declared and initialized X_BM...
const RevoluteJoint<double>& elbow =
plant.AddJoint<RevoluteJoint>(
"Elbow", /* joint name */
plant.world_body(), /* parent body */
{}, /* frame F IS the world frame W */
pendulum, /* child body, the pendulum */
X_BM, /* pose of frame M in the body frame B */
Vector3d::UnitZ()); /* revolute axis in this case */
Warning
Do not ever attempt to instantiate and manipulate Joint objects on the stack; it will fail. Add joints to your plant using the provided API MultibodyPlant::AddJoint() as in the example above.
Template Parameters
TThe scalar type, which must be one of the default scalars.

#include <drake/multibody/tree/multibody_tree.h>

Classes

struct  BluePrint
 (Advanced) Structure containing all the information needed to build the MultibodyTree implementation for a Joint. More...
 
struct  JointImplementation
 (Advanced) A Joint is implemented in terms of MultibodyTree elements, typically a Mobilizer. More...
 

Public Member Functions

 Joint (const std::string &name, const Frame< T > &frame_on_parent, const Frame< T > &frame_on_child, VectorX< double > damping, const VectorX< double > &pos_lower_limits, const VectorX< double > &pos_upper_limits, const VectorX< double > &vel_lower_limits, const VectorX< double > &vel_upper_limits, const VectorX< double > &acc_lower_limits, const VectorX< double > &acc_upper_limits)
 Creates a joint between two Frame objects which imposes a given kinematic relation between frame F attached on the parent body P and frame M attached on the child body B. More...
 
 Joint (const std::string &name, const Frame< T > &frame_on_parent, const Frame< T > &frame_on_child, const VectorX< double > &pos_lower_limits, const VectorX< double > &pos_upper_limits, const VectorX< double > &vel_lower_limits, const VectorX< double > &vel_upper_limits, const VectorX< double > &acc_lower_limits, const VectorX< double > &acc_upper_limits)
 Additional constructor overload for joints with zero damping. More...
 
virtual ~Joint ()
 
JointIndex index () const
 Returns this element's unique index. More...
 
const std::string & name () const
 Returns the name of this joint. More...
 
const RigidBody< T > & parent_body () const
 Returns a const reference to the parent body P. More...
 
const RigidBody< T > & child_body () const
 Returns a const reference to the child body B. More...
 
const Frame< T > & frame_on_parent () const
 Returns a const reference to the frame F attached on the parent body P. More...
 
const Frame< T > & frame_on_child () const
 Returns a const reference to the frame M attached on the child body B. More...
 
virtual const std::string & type_name () const =0
 Returns a string identifying the type of this joint, such as "revolute" or "prismatic". More...
 
int velocity_start () const
 Returns the index to the first generalized velocity for this joint within the vector v of generalized velocities for the full multibody system. More...
 
int num_velocities () const
 Returns the number of generalized velocities describing this joint. More...
 
int position_start () const
 Returns the index to the first generalized position for this joint within the vector q of generalized positions for the full multibody system. More...
 
int num_positions () const
 Returns the number of generalized positions describing this joint. More...
 
bool can_rotate () const
 Returns true if this joint's mobility allows relative rotation of the two frames associated with this joint. More...
 
bool can_translate () const
 Returns true if this joint's mobility allows relative translation of the two frames associated with this joint. More...
 
std::string position_suffix (int position_index_in_joint) const
 Returns a string suffix (e.g. More...
 
std::string velocity_suffix (int velocity_index_in_joint) const
 Returns a string suffix (e.g. More...
 
const T & GetOnePosition (const systems::Context< T > &context) const
 Returns the position coordinate for joints with a single degree of freedom. More...
 
const T & GetOneVelocity (const systems::Context< T > &context) const
 Returns the velocity coordinate for joints with a single degree of freedom. More...
 
void AddInOneForce (const systems::Context< T > &context, int joint_dof, const T &joint_tau, MultibodyForces< T > *forces) const
 Adds into forces a force along the one of the joint's degrees of freedom indicated by index joint_dof. More...
 
void AddInDamping (const systems::Context< T > &context, MultibodyForces< T > *forces) const
 Adds into forces the force due to damping within this joint. More...
 
void Lock (systems::Context< T > *context) const
 Lock the joint. More...
 
void Unlock (systems::Context< T > *context) const
 Unlock the joint. More...
 
bool is_locked (const systems::Context< T > &context) const
 
const VectorX< double > & default_damping_vector () const
 Returns all default damping coefficients for joints that model viscous damping, of size num_velocities(). More...
 
const VectorX< double > & damping_vector () const
 (Deprecated.) More...
 
const VectorX< T > & GetDampingVector (const systems::Context< T > &context) const
 Returns the Context dependent damping coefficients stored as parameters in context. More...
 
void set_default_damping_vector (const VectorX< double > &damping)
 Sets the default value of the viscous damping coefficients for this joint. More...
 
void SetDampingVector (systems::Context< T > *context, const VectorX< T > &damping) const
 Sets the value of the viscous damping coefficients for this joint, stored as parameters in context. More...
 
Does not allow copy, move, or assignment
 Joint (const Joint &)=delete
 
Jointoperator= (const Joint &)=delete
 
 Joint (Joint &&)=delete
 
Jointoperator= (Joint &&)=delete
 
Methods to get and set the limits of <tt>this</tt> joint. For position

limits, the layout is the same as the generalized position's.

For velocity and acceleration limits, the layout is the same as the generalized velocity's. A limit with value +/- ∞ implies no upper or lower limit.

Returns the position lower limits.

const VectorX< double > & position_lower_limits () const
 
const VectorX< double > & position_upper_limits () const
 Returns the position upper limits. More...
 
const VectorX< double > & velocity_lower_limits () const
 Returns the velocity lower limits. More...
 
const VectorX< double > & velocity_upper_limits () const
 Returns the velocity upper limits. More...
 
const VectorX< double > & acceleration_lower_limits () const
 Returns the acceleration lower limits. More...
 
const VectorX< double > & acceleration_upper_limits () const
 Returns the acceleration upper limits. More...
 
const VectorX< double > & default_positions () const
 Returns the default positions. More...
 
void set_position_limits (const VectorX< double > &lower_limits, const VectorX< double > &upper_limits)
 Sets the position limits to lower_limits and upper_limits. More...
 
void set_velocity_limits (const VectorX< double > &lower_limits, const VectorX< double > &upper_limits)
 Sets the velocity limits to lower_limits and upper_limits. More...
 
void set_acceleration_limits (const VectorX< double > &lower_limits, const VectorX< double > &upper_limits)
 Sets the acceleration limits to lower_limits and upper_limits. More...
 
void set_default_positions (const VectorX< double > &default_positions)
 Sets the default positions to default_positions. More...
 
void SetDefaultPose (const math::RigidTransform< double > &X_FM)
 Sets this Joint's default generalized positions q₀ such that the pose of the child frame M in the parent frame F best matches the given pose. More...
 
math::RigidTransform< doubleGetDefaultPose () const
 Returns this Joint's default pose as a RigidTransform X_FM. More...
 
void SetDefaultPosePair (const Quaternion< double > &q_FM, const Vector3< double > &p_FM)
 (Advanced) This is the same as SetDefaultPose() except it takes the pose as a (quaternion, translation vector) pair. More...
 
std::pair< Eigen::Quaternion< double >, Vector3< double > > GetDefaultPosePair () const
 (Advanced) This is the same as GetDefaultPose() except it returns this Joint's default pose as a (quaternion, translation vector) pair. More...
 
- Public Member Functions inherited from MultibodyElement< T >
virtual ~MultibodyElement ()
 
ModelInstanceIndex model_instance () const
 Returns the ModelInstanceIndex of the model instance to which this element belongs. More...
 
template<typename MultibodyPlantDeferred = MultibodyPlant<T>>
const MultibodyPlantDeferred & GetParentPlant () const
 Returns the MultibodyPlant that owns this MultibodyElement. More...
 
void DeclareParameters (internal::MultibodyTreeSystem< T > *tree_system)
 Declares MultibodyTreeSystem Parameters at MultibodyTreeSystem::Finalize() time. More...
 
void SetDefaultParameters (systems::Parameters< T > *parameters) const
 Sets default values of parameters belonging to each MultibodyElement in parameters at a call to MultibodyTreeSystem::SetDefaultParameters(). More...
 
 MultibodyElement (const MultibodyElement &)=delete
 
MultibodyElementoperator= (const MultibodyElement &)=delete
 
 MultibodyElement (MultibodyElement &&)=delete
 
MultibodyElementoperator= (MultibodyElement &&)=delete
 

Protected Member Functions

virtual int do_get_velocity_start () const =0
 Implementation of the NVI velocity_start(), see velocity_start() for details. More...
 
virtual int do_get_num_velocities () const =0
 Implementation of the NVI num_velocities(), see num_velocities() for details. More...
 
virtual int do_get_position_start () const =0
 Implementation of the NVI position_start(), see position_start() for details. More...
 
virtual int do_get_num_positions () const =0
 Implementation of the NVI num_positions(), see num_positions() for details. More...
 
virtual std::string do_get_position_suffix (int index) const =0
 Implementation of the NVI position_suffix(), see position_suffix() for details. More...
 
virtual std::string do_get_velocity_suffix (int index) const =0
 Implementation of the NVI velocity_suffix(), see velocity_suffix() for details. More...
 
virtual void do_set_default_positions (const VectorX< double > &default_positions)=0
 Implementation of the NVI set_default_positions(), see set_default_positions() for details. More...
 
virtual void DoSetDefaultPosePair (const Quaternion< double > &q_FM, const Vector3< double > &p_FM)
 Implementation of the NVI SetDefaultPose(). More...
 
virtual std::pair< Eigen::Quaternion< double >, Vector3< double > > DoGetDefaultPosePair () const
 Implementation of the NVI GetDefaultPose(). More...
 
virtual const T & DoGetOnePosition (const systems::Context< T > &) const
 Implementation of the NVI GetOnePosition() that must only be implemented by those joint subclasses that have a single degree of freedom. More...
 
virtual const T & DoGetOneVelocity (const systems::Context< T > &) const
 Implementation of the NVI GetOneVelocity() that must only be implemented by those joint subclasses that have a single degree of freedom. More...
 
virtual void DoAddInOneForce (const systems::Context< T > &context, int joint_dof, const T &joint_tau, MultibodyForces< T > *forces) const =0
 Adds into forces a force along the one of the joint's degrees of freedom given by joint_dof. More...
 
virtual void DoAddInDamping (const systems::Context< T > &, MultibodyForces< T > *) const
 Adds into MultibodyForces the forces due to damping within this joint. More...
 
void DoSetTopology (const internal::MultibodyTreeTopology &) override
 Implementation of the NVI SetTopology(). More...
 
virtual std::unique_ptr< BluePrintMakeImplementationBlueprint () const =0
 This method must be implemented by derived classes in order to provide JointImplementationBuilder a BluePrint of their internal implementation JointImplementation. More...
 
const JointImplementationget_implementation () const
 Returns a const reference to the internal implementation of this joint. More...
 
bool has_implementation () const
 Returns whether this joint owns a particular implementation. More...
 
Methods to make a clone templated on different scalar types.

Clones this Joint (templated on T) to a joint templated on double.

virtual std::unique_ptr< Joint< double > > DoCloneToScalar (const internal::MultibodyTree< double > &tree_clone) const =0
 
virtual std::unique_ptr< Joint< AutoDiffXd > > DoCloneToScalar (const internal::MultibodyTree< AutoDiffXd > &tree_clone) const =0
 Clones this Joint (templated on T) to a joint templated on AutoDiffXd. More...
 
virtual std::unique_ptr< Joint< symbolic::Expression > > DoCloneToScalar (const internal::MultibodyTree< symbolic::Expression > &) const =0
 
- Protected Member Functions inherited from MultibodyElement< T >
 MultibodyElement ()
 Default constructor made protected so that sub-classes can still declare their default constructors if they need to. More...
 
 MultibodyElement (ModelInstanceIndex model_instance)
 Constructor which allows specifying a model instance. More...
 
 MultibodyElement (ModelInstanceIndex model_instance, int64_t index)
 Both the model instance and element index are specified. More...
 
template<typename ElementIndexType >
ElementIndexType index_impl () const
 Returns this element's unique index. More...
 
const internal::MultibodyTree< T > & get_parent_tree () const
 Returns a constant reference to the parent MultibodyTree that owns this element. More...
 
const internal::MultibodyTreeSystem< T > & GetParentTreeSystem () const
 Returns a constant reference to the parent MultibodyTreeSystem that owns the parent MultibodyTree that owns this element. More...
 
void SetTopology (const internal::MultibodyTreeTopology &tree)
 Gives MultibodyElement-derived objects the opportunity to retrieve their topology after MultibodyTree::Finalize() is invoked. More...
 
systems::NumericParameterIndex DeclareNumericParameter (internal::MultibodyTreeSystem< T > *tree_system, const systems::BasicVector< T > &model_vector)
 To be used by MultibodyElement-derived objects when declaring parameters in their implementation of DoDeclareParameters(). More...
 
systems::AbstractParameterIndex DeclareAbstractParameter (internal::MultibodyTreeSystem< T > *tree_system, const AbstractValue &model_value)
 To be used by MultibodyElement-derived objects when declaring parameters in their implementation of DoDeclareParameters(). More...
 

Friends

template<typename >
class Joint
 
class internal::JointImplementationBuilder< T >
 

Constructor & Destructor Documentation

◆ Joint() [1/4]

Joint ( const Joint< T > &  )
delete

◆ Joint() [2/4]

Joint ( Joint< T > &&  )
delete

◆ Joint() [3/4]

Joint ( const std::string &  name,
const Frame< T > &  frame_on_parent,
const Frame< T > &  frame_on_child,
VectorX< double damping,
const VectorX< double > &  pos_lower_limits,
const VectorX< double > &  pos_upper_limits,
const VectorX< double > &  vel_lower_limits,
const VectorX< double > &  vel_upper_limits,
const VectorX< double > &  acc_lower_limits,
const VectorX< double > &  acc_upper_limits 
)

Creates a joint between two Frame objects which imposes a given kinematic relation between frame F attached on the parent body P and frame M attached on the child body B.

The joint will be initialized to the model instance from frame_on_child (this is the typical convention for joints between the world and a model, or between two models (e.g. an arm to a gripper)). See this class's documentation for further details.

Parameters
[in]nameA string with a name identifying this joint.
[in]frame_on_parentThe frame F attached on the parent body connected by this joint.
[in]frame_on_childThe frame M attached on the child body connected by this joint.
[in]dampingA vector of viscous damping coefficients, of size num_velocities(). See default_damping_vector() for details.
[in]pos_lower_limitsA vector storing the lower limit for each generalized position. It must have the same size as pos_upper_limit. A value equal to -∞ implies no lower limit.
[in]pos_upper_limitsA vector storing the upper limit for each generalized position. It must have the same size as pos_lower_limit. A value equal to +∞ implies no upper limit.
[in]vel_lower_limitsA vector storing the lower limit for each generalized velocity. It must have the same size as vel_upper_limit. A value equal to -∞ implies no lower limit.
[in]vel_upper_limitsA vector storing the upper limit for each generalized velocity. It must have the same size as vel_lower_limit. A value equal to +∞ implies no upper limit.
[in]acc_lower_limitsA vector storing the lower limit for each generalized acceleration. It must have the same size as acc_upper_limit. A value equal to -∞ implies no lower limit.
[in]acc_upper_limitsA vector storing the upper limit for each generalized acceleration. It must have the same size as acc_lower_limit. A value equal to +∞ implies no upper limit.

◆ Joint() [4/4]

Joint ( const std::string &  name,
const Frame< T > &  frame_on_parent,
const Frame< T > &  frame_on_child,
const VectorX< double > &  pos_lower_limits,
const VectorX< double > &  pos_upper_limits,
const VectorX< double > &  vel_lower_limits,
const VectorX< double > &  vel_upper_limits,
const VectorX< double > &  acc_lower_limits,
const VectorX< double > &  acc_upper_limits 
)

Additional constructor overload for joints with zero damping.

Refer to the more general constructor signature taking damping for further details on the rest of the arguments for this constructor.

◆ ~Joint()

virtual ~Joint ( )
virtual

Member Function Documentation

◆ acceleration_lower_limits()

const VectorX<double>& acceleration_lower_limits ( ) const

Returns the acceleration lower limits.

◆ acceleration_upper_limits()

const VectorX<double>& acceleration_upper_limits ( ) const

Returns the acceleration upper limits.

◆ AddInDamping()

void AddInDamping ( const systems::Context< T > &  context,
MultibodyForces< T > *  forces 
) const

Adds into forces the force due to damping within this joint.

Parameters
[in]contextThe context storing the state and parameters for the model to which this joint belongs.
[out]forcesOn return, this method will add the force due to damping within this joint. This method aborts if forces is nullptr or if forces does not have the right sizes to accommodate a set of forces for the model to which this joint belongs.

◆ AddInOneForce()

void AddInOneForce ( const systems::Context< T > &  context,
int  joint_dof,
const T &  joint_tau,
MultibodyForces< T > *  forces 
) const

Adds into forces a force along the one of the joint's degrees of freedom indicated by index joint_dof.

The meaning for this degree of freedom and even its dimensional units depend on the specific joint sub-class. For a RevoluteJoint for instance, joint_dof can only be 0 since revolute joints's motion subspace only has one degree of freedom, while the units of joint_tau are those of torque (N⋅m in the MKS system of units). For multi-dof joints please refer to the documentation provided by specific joint sub-classes regarding the meaning of joint_dof.

Parameters
[in]contextThe context storing the state and parameters for the model to which this joint belongs.
[in]joint_dofIndex specifying one of the degrees of freedom for this joint. The index must be in the range 0 <= joint_dof < num_velocities() or otherwise this method will abort.
[in]joint_tauGeneralized force corresponding to the degree of freedom indicated by joint_dof to be added into forces.
[out]forcesOn return, this method will add force joint_tau for the degree of freedom joint_dof into the output forces. This method aborts if forces is nullptr or if forces doest not have the right sizes to accommodate a set of forces for the model to which this joint belongs.

◆ can_rotate()

bool can_rotate ( ) const

Returns true if this joint's mobility allows relative rotation of the two frames associated with this joint.

Precondition
the MultibodyPlant must be finalized.
See also
can_translate()

◆ can_translate()

bool can_translate ( ) const

Returns true if this joint's mobility allows relative translation of the two frames associated with this joint.

Precondition
the MultibodyPlant must be finalized.
See also
can_rotate()

◆ child_body()

const RigidBody<T>& child_body ( ) const

Returns a const reference to the child body B.

◆ damping_vector()

const VectorX<double>& damping_vector ( ) const

(Deprecated.)

Deprecated:
Use default_damping_vector() instead.
This will be removed from Drake on or after 2024-06-01.

◆ default_damping_vector()

const VectorX<double>& default_damping_vector ( ) const

Returns all default damping coefficients for joints that model viscous damping, of size num_velocities().

Joints that do not model damping return a zero vector of size num_velocities(). If vj is the vector of generalized velocities for this joint, of size num_velocities(), viscous damping models a generalized force at the joint of the form tau = -diag(dj)⋅vj, with dj the vector returned by this function. The units of the coefficients will depend on the specific joint type. For instance, for a revolute joint where vj is an angular velocity with units of rad/s and tau having units of N⋅m, the coefficient of viscous damping has units of N⋅m⋅s. Refer to each joint's documentation for further details.

◆ default_positions()

const VectorX<double>& default_positions ( ) const

Returns the default positions.

◆ do_get_num_positions()

virtual int do_get_num_positions ( ) const
protectedpure virtual

Implementation of the NVI num_positions(), see num_positions() for details.

Note
Implementations must meet the styleguide requirements for snake_case accessor methods.

◆ do_get_num_velocities()

virtual int do_get_num_velocities ( ) const
protectedpure virtual

Implementation of the NVI num_velocities(), see num_velocities() for details.

Note
Implementations must meet the styleguide requirements for snake_case accessor methods.

◆ do_get_position_start()

virtual int do_get_position_start ( ) const
protectedpure virtual

Implementation of the NVI position_start(), see position_start() for details.

Note
Implementations must meet the styleguide requirements for snake_case accessor methods.

◆ do_get_position_suffix()

virtual std::string do_get_position_suffix ( int  index) const
protectedpure virtual

Implementation of the NVI position_suffix(), see position_suffix() for details.

The suffix should contain only alphanumeric characters (e.g. 'wx' not '_wx' or '.wx').

◆ do_get_velocity_start()

virtual int do_get_velocity_start ( ) const
protectedpure virtual

Implementation of the NVI velocity_start(), see velocity_start() for details.

Note that this must be the offset within just the velocity vector, not within the composite state vector.

Note
Implementations must meet the styleguide requirements for snake_case accessor methods.

◆ do_get_velocity_suffix()

virtual std::string do_get_velocity_suffix ( int  index) const
protectedpure virtual

Implementation of the NVI velocity_suffix(), see velocity_suffix() for details.

The suffix should contain only alphanumeric characters (e.g. 'wx' not '_wx' or '.wx').

◆ do_set_default_positions()

virtual void do_set_default_positions ( const VectorX< double > &  default_positions)
protectedpure virtual

Implementation of the NVI set_default_positions(), see set_default_positions() for details.

It is the responsibility of the subclass to ensure that its joint implementation (i.e., mobilizer), should it have one, is updated with default_positions. Note that the Joint base class also stores default_positions (as a VectorX); the implementing mobilizer should have the same value but as a fixed-size vector.

Note
Implementations must meet the styleguide requirements for snake_case accessor methods.

◆ DoAddInDamping()

virtual void DoAddInDamping ( const systems::Context< T > &  ,
MultibodyForces< T > *   
) const
protectedvirtual

Adds into MultibodyForces the forces due to damping within this joint.

How forces are added to a MultibodyTree model depends on the underlying implementation of a particular joint (for instance, mobilizer vs. constraint) and therefore specific Joint subclasses must provide a definition for this method. The default implementation is a no-op for joints with no damping.

Reimplemented in QuaternionFloatingJoint< T >, RpyFloatingJoint< T >, RevoluteJoint< T >, PrismaticJoint< T >, BallRpyJoint< T >, and UniversalJoint< T >.

◆ DoAddInOneForce()

virtual void DoAddInOneForce ( const systems::Context< T > &  context,
int  joint_dof,
const T &  joint_tau,
MultibodyForces< T > *  forces 
) const
protectedpure virtual

Adds into forces a force along the one of the joint's degrees of freedom given by joint_dof.

How forces are added to a MultibodyTree model depends on the underlying implementation of a particular joint and therefore specific Joint subclasses must provide a definition for this method. For instance, a revolute joint could be modeled with a single generalized coordinate for the angular rotation (implemented through a RevoluteMobilizer) or it could be modeled using a constraint that only allows rotation about the joint's axis but that constrains the motion in the other five degrees of freedom. This method is only called by the public NVI AddInOneForce() and therefore input arguments were checked to be valid.

See also
The public NVI AddInOneForce() for details.

Implemented in QuaternionFloatingJoint< T >, RpyFloatingJoint< T >, RevoluteJoint< T >, PrismaticJoint< T >, BallRpyJoint< T >, UniversalJoint< T >, and WeldJoint< T >.

◆ DoCloneToScalar() [1/3]

virtual std::unique_ptr<Joint<double> > DoCloneToScalar ( const internal::MultibodyTree< double > &  tree_clone) const
protectedpure virtual

◆ DoCloneToScalar() [2/3]

virtual std::unique_ptr<Joint<AutoDiffXd> > DoCloneToScalar ( const internal::MultibodyTree< AutoDiffXd > &  tree_clone) const
protectedpure virtual

Clones this Joint (templated on T) to a joint templated on AutoDiffXd.

◆ DoCloneToScalar() [3/3]

virtual std::unique_ptr<Joint<symbolic::Expression> > DoCloneToScalar ( const internal::MultibodyTree< symbolic::Expression > &  ) const
protectedpure virtual

◆ DoGetDefaultPosePair()

virtual std::pair<Eigen::Quaternion<double>, Vector3<double> > DoGetDefaultPosePair ( ) const
protectedvirtual

Implementation of the NVI GetDefaultPose().

This is optional for Joint subclasses except for floating (6 dof) Joints. The subclass should convert its default_positions to pose X_FM and return that as a (quaternion, translation) pair. If the subclass already uses (quaternion, translation) as generalized coordinates (i.e. it's a quaternion_floating_joint) it must return those exactly (don't convert to a transform first).

◆ DoGetOnePosition()

virtual const T& DoGetOnePosition ( const systems::Context< T > &  ) const
protectedvirtual

Implementation of the NVI GetOnePosition() that must only be implemented by those joint subclasses that have a single degree of freedom.

The default implementation for all other joints is to abort with an appropriate message. Revolute and prismatic are examples of joints that will want to implement this method.

◆ DoGetOneVelocity()

virtual const T& DoGetOneVelocity ( const systems::Context< T > &  ) const
protectedvirtual

Implementation of the NVI GetOneVelocity() that must only be implemented by those joint subclasses that have a single degree of freedom.

The default implementation for all other joints is to abort with an appropriate message. Revolute and prismatic are examples of joints that will want to implement this method.

◆ DoSetDefaultPosePair()

virtual void DoSetDefaultPosePair ( const Quaternion< double > &  q_FM,
const Vector3< double > &  p_FM 
)
protectedvirtual

Implementation of the NVI SetDefaultPose().

This is optional for Joint subclasses except for floating (6 dof) Joints. The subclass should convert the input to the closest equivalent in generalized coordinates and invoke set_default_positions() to record them. If the subclass already uses (quaternion, translation) as generalized coordinates (i.e. it's a quaternion_floating_joint) it must store those exactly.

◆ DoSetTopology()

void DoSetTopology ( const internal::MultibodyTreeTopology &  tree)
overrideprotectedvirtual

Implementation of the NVI SetTopology().

For advanced use only for developers implementing new MultibodyTree components.

Implements MultibodyElement< T >.

◆ frame_on_child()

const Frame<T>& frame_on_child ( ) const

Returns a const reference to the frame M attached on the child body B.

◆ frame_on_parent()

const Frame<T>& frame_on_parent ( ) const

Returns a const reference to the frame F attached on the parent body P.

◆ get_implementation()

const JointImplementation& get_implementation ( ) const
protected

Returns a const reference to the internal implementation of this joint.

Warning
The MultibodyTree model must have already been finalized, or this method will abort.

◆ GetDampingVector()

const VectorX<T>& GetDampingVector ( const systems::Context< T > &  context) const

Returns the Context dependent damping coefficients stored as parameters in context.

Refer to default_damping_vector() for details.

Parameters
[in]contextThe context storing the state and parameters for the model to which this joint belongs.

◆ GetDefaultPose()

math::RigidTransform<double> GetDefaultPose ( ) const

Returns this Joint's default pose as a RigidTransform X_FM.

Note
Currently this is implemented only for floating (6 dof) joints which can represent any pose.
Exceptions
std::exceptionif called for any Joint type that does not implement this function.
Return values
X_FMThe default pose as a rigid transform.
See also
get_default_positions() to see the generalized positions q₀ that this joint used to generate the returned transform.

◆ GetDefaultPosePair()

std::pair<Eigen::Quaternion<double>, Vector3<double> > GetDefaultPosePair ( ) const

(Advanced) This is the same as GetDefaultPose() except it returns this Joint's default pose as a (quaternion, translation vector) pair.

Return values
q_FM,p_FMThe default pose as a (quaternion, translation) pair.
See also
GetDefaultPose() for more information

◆ GetOnePosition()

const T& GetOnePosition ( const systems::Context< T > &  context) const

Returns the position coordinate for joints with a single degree of freedom.

Exceptions
std::exceptionif the joint does not have a single degree of freedom.

◆ GetOneVelocity()

const T& GetOneVelocity ( const systems::Context< T > &  context) const

Returns the velocity coordinate for joints with a single degree of freedom.

Exceptions
std::exceptionif the joint does not have a single degree of freedom.

◆ has_implementation()

bool has_implementation ( ) const
protected

Returns whether this joint owns a particular implementation.

If the MultibodyTree has been finalized, this will return true.

◆ index()

JointIndex index ( ) const

Returns this element's unique index.

◆ is_locked()

bool is_locked ( const systems::Context< T > &  context) const
Returns
true if the joint is locked, false otherwise.

◆ Lock()

void Lock ( systems::Context< T > *  context) const

Lock the joint.

Its generalized velocities will be 0 until it is unlocked.

◆ MakeImplementationBlueprint()

virtual std::unique_ptr<BluePrint> MakeImplementationBlueprint ( ) const
protectedpure virtual

This method must be implemented by derived classes in order to provide JointImplementationBuilder a BluePrint of their internal implementation JointImplementation.

◆ name()

const std::string& name ( ) const

Returns the name of this joint.

◆ num_positions()

int num_positions ( ) const

Returns the number of generalized positions describing this joint.

◆ num_velocities()

int num_velocities ( ) const

Returns the number of generalized velocities describing this joint.

◆ operator=() [1/2]

Joint& operator= ( Joint< T > &&  )
delete

◆ operator=() [2/2]

Joint& operator= ( const Joint< T > &  )
delete

◆ parent_body()

const RigidBody<T>& parent_body ( ) const

Returns a const reference to the parent body P.

◆ position_lower_limits()

const VectorX<double>& position_lower_limits ( ) const

◆ position_start()

int position_start ( ) const

Returns the index to the first generalized position for this joint within the vector q of generalized positions for the full multibody system.

◆ position_suffix()

std::string position_suffix ( int  position_index_in_joint) const

Returns a string suffix (e.g.

to be appended to the name()) to identify the kth position in this joint. position_index_in_joint must be in [0, num_positions()).

Precondition
the MultibodyPlant must be finalized.

◆ position_upper_limits()

const VectorX<double>& position_upper_limits ( ) const

Returns the position upper limits.

◆ set_acceleration_limits()

void set_acceleration_limits ( const VectorX< double > &  lower_limits,
const VectorX< double > &  upper_limits 
)

Sets the acceleration limits to lower_limits and upper_limits.

Exceptions
std::exceptionif the dimension of lower_limits or upper_limits does not match num_velocities().
std::exceptionif any of lower_limits is larger than the corresponding term in upper_limits.

◆ set_default_damping_vector()

void set_default_damping_vector ( const VectorX< double > &  damping)

Sets the default value of the viscous damping coefficients for this joint.

Refer to default_damping_vector() for details.

Exceptions
std::exceptionif damping.size() != num_velocities().
std::exceptionif any of the damping coefficients is negative.
Precondition
the MultibodyPlant must not be finalized.

◆ set_default_positions()

void set_default_positions ( const VectorX< double > &  default_positions)

Sets the default positions to default_positions.

Joint subclasses are expected to implement the do_set_default_positions().

Exceptions
std::exceptionif the dimension of default_positions does not match num_positions().
Note
The values in default_positions are NOT constrained to be within position_lower_limits() and position_upper_limits().

◆ set_position_limits()

void set_position_limits ( const VectorX< double > &  lower_limits,
const VectorX< double > &  upper_limits 
)

Sets the position limits to lower_limits and upper_limits.

Exceptions
std::exceptionif the dimension of lower_limits or upper_limits does not match num_positions().
std::exceptionif any of lower_limits is larger than the corresponding term in upper_limits.
Note
Setting the position limits does not affect the default_positions(), regardless of whether the current default_positions() satisfy the new position limits.

◆ set_velocity_limits()

void set_velocity_limits ( const VectorX< double > &  lower_limits,
const VectorX< double > &  upper_limits 
)

Sets the velocity limits to lower_limits and upper_limits.

Exceptions
std::exceptionif the dimension of lower_limits or upper_limits does not match num_velocities().
std::exceptionif any of lower_limits is larger than the corresponding term in upper_limits.

◆ SetDampingVector()

void SetDampingVector ( systems::Context< T > *  context,
const VectorX< T > &  damping 
) const

Sets the value of the viscous damping coefficients for this joint, stored as parameters in context.

Refer to default_damping_vector() for details.

Parameters
[out]contextThe context storing the state and parameters for the model to which this joint belongs.
[in]dampingThe vector of damping values.
Exceptions
std::exceptionif damping.size() != num_velocities().
std::exceptionif any of the damping coefficients is negative.
Note
Some multi-dof joints may have specific semantics for their damping vector that are not enforced here. For instance, QuaternionFloatingJoint assumes identical damping values for all 3 angular velocity components and identical damping values for all 3 translational velocity components. It will thus use angular_damping = damping[0] and translational_damping = damping[3]. Refer to the particular subclass for more semantic information.

◆ SetDefaultPose()

void SetDefaultPose ( const math::RigidTransform< double > &  X_FM)

Sets this Joint's default generalized positions q₀ such that the pose of the child frame M in the parent frame F best matches the given pose.

The pose is given by a RigidTransform X_FM, but a Joint will represent pose differently.

Note
Currently this is implemented only for floating (6 dof) joints which can represent any pose.
Exceptions
std::exceptionif called for any Joint type that does not implement this function.
See also
get_default_positions() to see the resulting q₀ after this call.
SetDefaultPosePair() for an alternative using a quaternion

◆ SetDefaultPosePair()

void SetDefaultPosePair ( const Quaternion< double > &  q_FM,
const Vector3< double > &  p_FM 
)

(Advanced) This is the same as SetDefaultPose() except it takes the pose as a (quaternion, translation vector) pair.

See also
SetDefaultPose() for more information

◆ type_name()

virtual const std::string& type_name ( ) const
pure virtual

Returns a string identifying the type of this joint, such as "revolute" or "prismatic".

Implemented in ScrewJoint< T >, RevoluteJoint< T >, QuaternionFloatingJoint< T >, RpyFloatingJoint< T >, UniversalJoint< T >, PlanarJoint< T >, BallRpyJoint< T >, PrismaticJoint< T >, and WeldJoint< T >.

◆ Unlock()

void Unlock ( systems::Context< T > *  context) const

Unlock the joint.

◆ velocity_lower_limits()

const VectorX<double>& velocity_lower_limits ( ) const

Returns the velocity lower limits.

◆ velocity_start()

int velocity_start ( ) const

Returns the index to the first generalized velocity for this joint within the vector v of generalized velocities for the full multibody system.

◆ velocity_suffix()

std::string velocity_suffix ( int  velocity_index_in_joint) const

Returns a string suffix (e.g.

to be appended to the name()) to identify the kth velocity in this joint. velocity_index_in_joint must be in [0, num_velocities()).

Precondition
the MultibodyPlant must be finalized.

◆ velocity_upper_limits()

const VectorX<double>& velocity_upper_limits ( ) const

Returns the velocity upper limits.

Friends And Related Function Documentation

◆ internal::JointImplementationBuilder< T >

friend class internal::JointImplementationBuilder< T >
friend

◆ Joint

friend class Joint
friend

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