Drake
MultibodyPlant< T > Class Template Reference

MultibodyPlant is a Drake system framework representation (see systems::System) for the model of a physical system consisting of a collection of interconnected bodies. More...

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

Public Member Functions

 MultibodyPlant (double time_step=0)
 Default constructor creates a plant with a single "world" body. More...
 
template<typename U >
 MultibodyPlant (const MultibodyPlant< U > &other)
 Scalar-converting copy constructor. See System Scalar Conversion. More...
 
int num_bodies () const
 Returns the number of bodies in the model, including the "world" body, which is always part of the model. More...
 
int num_joints () const
 Returns the number of joints in the model. More...
 
int num_actuators () const
 Returns the number of joint actuators in the model. More...
 
int num_model_instances () const
 Returns the number of model instances in the model. More...
 
int num_positions () const
 Returns the size of the generalized position vector q for this MultibodyPlant. More...
 
int num_positions (ModelInstanceIndex model_instance) const
 Returns the size of the generalized position vector q for a specific model instance. More...
 
int num_velocities () const
 Returns the size of the generalized velocity vector v for this MultibodyPlant. More...
 
int num_velocities (ModelInstanceIndex model_instance) const
 Returns the size of the generalized velocity vector v for a specific model instance. More...
 
int num_multibody_states () const
 Returns the size of the multibody system state vector x = [q; v] for this MultibodyPlant. More...
 
int num_actuated_dofs () const
 Returns the total number of actuated degrees of freedom. More...
 
int num_actuated_dofs (ModelInstanceIndex model_instance) const
 Returns the total number of actuated degrees of freedom for a specific model instance. More...
 
geometry::SourceId RegisterAsSourceForSceneGraph (geometry::SceneGraph< T > *scene_graph)
 Registers this plant to serve as a source for an instance of SceneGraph. More...
 
geometry::GeometryId RegisterVisualGeometry (const Body< T > &body, const Isometry3< double > &X_BG, const geometry::Shape &shape, const std::string &name, const geometry::VisualMaterial &material, geometry::SceneGraph< T > *scene_graph)
 Registers geometry in a SceneGraph with a given geometry::Shape to be used for visualization of a given body. More...
 
geometry::GeometryId RegisterVisualGeometry (const Body< T > &body, const Isometry3< double > &X_BG, const geometry::Shape &shape, const std::string &name, geometry::SceneGraph< T > *scene_graph)
 Overload for visual geometry registration; it implicitly assigns the default material. More...
 
const std::vector< geometry::GeometryId > & GetVisualGeometriesForBody (const Body< T > &body) const
 Returns an array of GeometryId's identifying the different visual geometries for body previously registered with a SceneGraph. More...
 
int num_visual_geometries () const
 Returns the number of geometries registered for visualization. More...
 
geometry::GeometryId RegisterCollisionGeometry (const Body< T > &body, const Isometry3< double > &X_BG, const geometry::Shape &shape, const std::string &name, const CoulombFriction< double > &coulomb_friction, geometry::SceneGraph< T > *scene_graph)
 Registers geometry in a SceneGraph with a given geometry::Shape to be used for the contact modeling of a given body. More...
 
const std::vector< geometry::GeometryId > & GetCollisionGeometriesForBody (const Body< T > &body) const
 Returns an array of GeometryId's identifying the different contact geometries for body previously registered with a SceneGraph. More...
 
int num_collision_geometries () const
 Returns the number of geometries registered for contact modeling. More...
 
geometry::GeometrySet CollectRegisteredGeometries (const std::vector< const Body< T > *> &bodies) const
 For each of the provided bodies, collects up all geometries that have been registered to that body. More...
 
const CoulombFriction< double > & default_coulomb_friction (geometry::GeometryId id) const
 Returns the friction coefficients provided during geometry registration for the given geometry id. More...
 
bool geometry_source_is_registered () const
 Returns true if this MultibodyPlant was registered with a SceneGraph. More...
 
const Body< T > * GetBodyFromFrameId (geometry::FrameId frame_id) const
 Given a geometry frame identifier, returns a pointer to the body associated with that id (nullptr if there is no such body). More...
 
optional< geometry::FrameIdGetBodyFrameIdIfExists (BodyIndex body_index) const
 If the body with body_index has geometry registered with it, it returns the geometry::FrameId associated with it. More...
 
geometry::FrameId GetBodyFrameIdOrThrow (BodyIndex body_index) const
 If the body with body_index has geometry registered with it, it returns the geometry::FrameId associated with it. More...
 
const systems::OutputPort< T > & get_generalized_contact_forces_output_port (ModelInstanceIndex model_instance) const
 Returns a constant reference to the output port of generalized contact forces for a specific model instance. More...
 
const systems::OutputPort< T > & get_contact_results_output_port () const
 Returns a constant reference to the port that outputs ContactResults. More...
 
const RigidBody< T > & world_body () const
 Returns a constant reference to the world body. More...
 
const BodyFrame< T > & world_frame () const
 Returns a constant reference to the world frame. More...
 
const MultibodyTree< T > & model () const
 Returns a constant reference to the underlying MultibodyTree model for this plant. More...
 
bool is_finalized () const
 Returns true if this MultibodyPlant was finalized with a call to Finalize(). More...
 
void Finalize (geometry::SceneGraph< T > *scene_graph=nullptr)
 This method must be called after all elements in the model (joints, bodies, force elements, constraints, etc.) are added and before any computations are performed. More...
 
double time_step () const
 The time step (or period) used to model this plant as a discrete system with periodic updates. More...
 
void SetDefaultState (const systems::Context< T > &context, systems::State< T > *state) const override
 Sets the state in context so that generalized positions and velocities are zero. More...
 
Does not allow copy, move, or assignment
 MultibodyPlant (const MultibodyPlant &)=delete
 
MultibodyPlantoperator= (const MultibodyPlant &)=delete
 
 MultibodyPlant (MultibodyPlant &&)=delete
 
MultibodyPlantoperator= (MultibodyPlant &&)=delete
 
Adding new multibody elements

MultibodyPlant users will add modeling elements like bodies, joints, force elements, constraints, etc, using one of these methods.

Once a user is done adding all modeling elements, the Finalize() method must be called before invoking any MultibodyPlant service to perform computations. An attempt to call any of these methods after a call to Finalize() on the plant, will result on an exception being thrown. See Finalize() for details.

const RigidBody< T > & AddRigidBody (const std::string &name, ModelInstanceIndex model_instance, const SpatialInertia< double > &M_BBo_B)
 Creates a rigid body with the provided name and spatial inertia. More...
 
const RigidBody< T > & AddRigidBody (const std::string &name, const SpatialInertia< double > &M_BBo_B)
 Creates a rigid body with the provided name and spatial inertia. More...
 
template<template< typename > class FrameType>
const FrameType< T > & AddFrame (std::unique_ptr< FrameType< T >> frame)
 This method adds a Frame of type FrameType<T>. More...
 
template<template< typename Scalar > class JointType>
const JointType< T > & AddJoint (std::unique_ptr< JointType< T >> joint)
 This method adds a Joint of type JointType between two bodies. More...
 
template<template< typename > class JointType, typename... Args>
const JointType< T > & AddJoint (const std::string &name, const Body< T > &parent, const optional< Isometry3< double >> &X_PF, const Body< T > &child, const optional< Isometry3< double >> &X_BM, Args &&... args)
 This method adds a Joint of type JointType between two bodies. More...
 
template<template< typename Scalar > class ForceElementType, typename... Args>
const ForceElementType< T > & AddForceElement (Args &&... args)
 Adds a new force element model of type ForceElementType to this MultibodyPlant. More...
 
template<template< typename Scalar > class ForceElementType, typename... Args>
std::enable_if< std::is_same< ForceElementType< T >, UniformGravityFieldElement< T > >::value, const ForceElementType< T > & >::type AddForceElement (Args &&... args)
 
const JointActuator< T > & AddJointActuator (const std::string &name, const Joint< T > &joint)
 Creates and adds a JointActuator model for an actuator acting on a given joint. More...
 
ModelInstanceIndex AddModelInstance (const std::string &name)
 Creates a new model instance. More...
 
const WeldJoint< T > & WeldFrames (const Frame< T > &A, const Frame< T > &B, const Isometry3< double > &X_AB=Isometry3< double >::Identity())
 Welds frames A and B with relative pose X_AB. More...
 
Querying for multibody elements by name

These methods allow a user to query whether a given multibody element is part of this plant's model.

These queries can be performed at any time during the lifetime of a MultibodyPlant model, i.e. there is no restriction on whether they must be called before or after Finalize(). That is, these queries can be performed while new multibody elements are being added to the model.

bool HasBodyNamed (const std::string &name) const
 
bool HasBodyNamed (const std::string &name, ModelInstanceIndex model_instance) const
 
bool HasJointNamed (const std::string &name) const
 
bool HasJointNamed (const std::string &name, ModelInstanceIndex model_instance) const
 
bool HasJointActuatorNamed (const std::string &name) const
 
bool HasJointActuatorNamed (const std::string &name, ModelInstanceIndex model_instance) const
 
bool HasModelInstanceNamed (const std::string &name) const
 
Retrieving multibody elements by name

These methods allow a user to retrieve a reference to a multibody element by its name.

An exception is thrown if there is no element with the requested name.

const Body< T > & GetBodyByName (const std::string &name) const
 These queries can be performed at any time during the lifetime of a MultibodyPlant, i.e. More...
 
const Body< T > & GetBodyByName (const std::string &name, ModelInstanceIndex model_instance) const
 Returns a constant reference to the body that is uniquely identified by the string name and model_instance in this MultibodyPlant. More...
 
const Frame< T > & GetFrameByName (const std::string &name) const
 Returns a constant reference to a frame that is identified by the string name in this model. More...
 
const Frame< T > & GetFrameByName (const std::string &name, ModelInstanceIndex model_instance) const
 Returns a constant reference to the frame that is uniquely identified by the string name in model_instance. More...
 
const Joint< T > & GetJointByName (const std::string &name) const
 Returns a constant reference to a joint that is identified by the string name in this MultibodyPlant. More...
 
const Joint< T > & GetJointByName (const std::string &name, ModelInstanceIndex model_instance) const
 Returns a constant reference to the joint that is uniquely identified by the string name and model_instance in this MultibodyPlant. More...
 
template<template< typename > class JointType>
const JointType< T > & GetJointByName (const std::string &name) const
 A templated version of GetJointByName() to return a constant reference of the specified type JointType in place of the base Joint class. More...
 
template<template< typename > class JointType>
const JointType< T > & GetJointByName (const std::string &name, ModelInstanceIndex model_instance) const
 A templated version of GetJointByName() to return a constant reference of the specified type JointType in place of the base Joint class. More...
 
const JointActuator< T > & GetJointActuatorByName (const std::string &name) const
 Returns a constant reference to an actuator that is identified by the string name in this MultibodyPlant. More...
 
const JointActuator< T > & GetJointActuatorByName (const std::string &name, ModelInstanceIndex model_instance) const
 Returns a constant reference to the actuator that is uniquely identified by the string name and model_instance in this MultibodyPlant. More...
 
ModelInstanceIndex GetModelInstanceByName (const std::string &name) const
 Returns the index to the model instance that is uniquely identified by the string name in this MultibodyPlant. More...
 
Retrieving ports for communication with a SceneGraph.
optional< geometry::SourceIdget_source_id () const
 Returns the unique id identifying this plant as a source for a SceneGraph. More...
 
const systems::InputPort< T > & get_geometry_query_input_port () const
 Returns a constant reference to the input port used to perform geometric queries on a SceneGraph. More...
 
const systems::OutputPort< T > & get_geometry_poses_output_port () const
 Returns the output port of frames' poses to communicate with a SceneGraph. More...
 
Actuation input

The input vector of actuation values can be provided either as a single input port which describes the entire plant (in the case where only a single model instance has actuated dofs), or through multiple input ports which each provide the actuation values for a specific model instance.

See AddJointActuator() and num_actuators().

const systems::InputPort< T > & get_actuation_input_port () const
 Returns a constant reference to the input port for external actuation for the case where only one model instance has actuated dofs. More...
 
const systems::InputPort< T > & get_actuation_input_port (ModelInstanceIndex model_instance) const
 Returns a constant reference to the input port for external actuation for a specific model instance. More...
 
Continuous state output

Output ports are provided to access the continuous state of the whole plant and for individual model instances.

const systems::OutputPort< T > & get_continuous_state_output_port () const
 Returns a constant reference to the output port for the full continuous state of the model. More...
 
const systems::OutputPort< T > & get_continuous_state_output_port (ModelInstanceIndex model_instance) const
 Returns a constant reference to the output port for the continuous state of a specific model instance. More...
 
Contact by penalty method

Currently MultibodyPlant uses a rigid contact model that is, bodies in the model are infinitely stiff or ideal rigid bodies. Therefore, the mathematical description of the rigid contact model needs to include non-penetration constraints among bodies in the formulation. There are several numerical methods to impose and solve these constraints. In a penalty method approach, we allow for a certain amount of interpenetration and we compute contact forces according to a simple law of the form:

  fₙ = k(1+dẋ)x

where the normal contact force fₙ is made a continuous function of the penetration distance x between the bodies (defined to be positive when the bodies are in contact) and the penetration distance rate ẋ (with ẋ > 0 meaning the penetration distance is increasing and therefore the interpenetration between the bodies is also increasing). k and d are the penalty method coefficients for stiffness and damping. These are ad-hoc parameters which need to be tuned as a trade-off between:

  • The accuracy of the numerical approximation to rigid contact, which requires a stiffness that approaches infinity, and
  • the computational cost of the numerical integration, which will require smaller time steps for stiffer systems.

There is no exact procedure for choosing these coefficients, and estimating them manually can be cumbersome since in general they will depend on the scale of the problem including masses, speeds and even body sizes. However, MultibodyPlant aids the estimation of these coefficients using a heuristic function based on a user-supplied "penetration allowance", see set_penetration_allowance(). The penetration allowance is a number in meters that specifies the order of magnitude of the average penetration between bodies in the system that the user is willing to accept as reasonable for the problem being solved. For instance, in the robotics manipulation of ordinary daily objects the user might set this number to 1 millimeter. However, the user might want to increase it for the simulation of heavy walking robots for which an allowance of 1 millimeter would result in a very stiff system.

As for the damping coefficient in the simple law above, MultibodyPlant chooses the damping coefficient d to model inelastic collisions and therefore sets it so that the penetration distance x behaves as a critically damped oscillator. That is, at the limit of ideal rigid contact (very stiff penalty coefficient k or equivalently the penetration allowance goes to zero), this method behaves as a unilateral constraint on the penetration distance, which models a perfect inelastic collision. For most applications, such as manipulation and walking, this is the desired behavior.

When set_penetration_allowance() is called, MultibodyPlant will estimate reasonable penalty method coefficients as a function of the input penetration allowance. Users will want to run their simulation a number of times and asses they are satisfied with the level of inter-penetration actually observed in the simulation; if the observed penetration is too large, the user will want to set a smaller penetration allowance. If the system is too stiff and the time integration requires very small time steps while at the same time the user can afford larger inter-penetrations, the user will want to increase the penetration allowance. Typically, the observed penetration will be proportional to the penetration allowance. Thus scaling the penetration allowance by say a factor of 0.5, would typically results in inter-penetrations being reduced by the same factor of 0.5. In summary, users should choose the largest penetration allowance that results in inter-penetration levels that are acceptable for the particular application (even when in theory this penetration should be zero for perfectly rigid bodies.)

For a given penetration allowance, the contact interaction that takes two bodies with a non-zero approaching velocity to zero approaching velocity, takes place in a finite amount of time (for ideal rigid contact this time is zero.) A good estimate of this time period is given by a call to get_contact_penalty_method_time_scale(). Users might want to query this value to either set the maximum time step in error-controlled time integration or to set the time step for fixed time step integration. As a guidance, typical fixed time step integrators will become unstable for time steps larger than about a tenth of this time scale.

void set_penetration_allowance (double penetration_allowance=0.001)
 Sets the penetration allowance used to estimate the coefficients in the penalty method used to impose non-penetration among bodies. More...
 
double get_contact_penalty_method_time_scale () const
 Returns a time-scale estimate tc based on the requested penetration allowance δ set with set_penetration_allowance(). More...
 
Stribeck model of friction

Currently MultibodyPlant uses the Stribeck approximation to model dry friction. The Stribeck model of friction is an approximation to Coulomb's law of friction that allows using continuous time integration without the need to specify complementarity constraints. While this results in a simpler model immediately tractable with standard numerical methods for integration of ODE's, it often leads to stiff dynamics that require an explicit integrator to take very small time steps. It is therefore recommended to use error controlled integrators when using this model. See tangent_force for a detailed discussion of the Stribeck model.

void set_stiction_tolerance (double v_stiction=0.001)
 Sets the stiction tolerance v_stiction for the Stribeck model, where v_stiction must be specified in m/s (meters per second.) v_stiction defaults to a value of 1 millimeter per second. More...
 
- Public Member Functions inherited from MultibodyTreeSystem< T >
 MultibodyTreeSystem (std::unique_ptr< MultibodyTree< T >> tree, bool is_discrete=false)
 Takes ownership of the given tree, finalizes it if it hasn't already been finalized, and then allocates the resources it needs. More...
 
template<typename U >
 MultibodyTreeSystem (const MultibodyTreeSystem< U > &other)
 Scalar-converting copy constructor. More...
 
 ~MultibodyTreeSystem () override
 
bool is_discrete () const
 
const MultibodyTree< T > & tree () const
 Returns a const reference to the MultibodyTree owned by this class. More...
 
const PositionKinematicsCache< T > & EvalPositionKinematics (const systems::Context< T > &context) const
 Returns a reference to the up to date PositionKinematicsCache in the given Context, recalculating it first if necessary. More...
 
const VelocityKinematicsCache< T > & EvalVelocityKinematics (const systems::Context< T > &context) const
 Returns a reference to the up to date VelocityKinematicsCache in the given Context, recalculating it first if necessary. More...
 
 MultibodyTreeSystem (const MultibodyTreeSystem &)=delete
 
MultibodyTreeSystemoperator= (const MultibodyTreeSystem &)=delete
 
 MultibodyTreeSystem (MultibodyTreeSystem &&)=delete
 
MultibodyTreeSystemoperator= (MultibodyTreeSystem &&)=delete
 
- Public Member Functions inherited from LeafSystem< T >
 ~LeafSystem () override
 
std::unique_ptr< CompositeEventCollection< T > > AllocateCompositeEventCollection () const final
 Allocates a CompositeEventCollection object for this system. More...
 
std::unique_ptr< LeafContext< T > > AllocateContext () const
 Shadows System<T>::AllocateContext to provide a more concrete return type LeafContext<T>. More...
 
std::unique_ptr< ContextBaseDoAllocateContext () const final
 Derived class implementations should allocate a suitable concrete Context type, then invoke the above InitializeContextBase() method. More...
 
void SetDefaultState (const Context< T > &context, State< T > *state) const override
 Default implementation: sets all continuous state to the model vector given in DeclareContinuousState (or zero if no model vector was given) and discrete states to zero. More...
 
void SetDefaultParameters (const Context< T > &context, Parameters< T > *parameters) const override
 Default implementation: sets all numeric parameters to the model vector given to DeclareNumericParameter, or else if no model was provided sets the numeric parameter to one. More...
 
std::unique_ptr< ContinuousState< T > > AllocateTimeDerivatives () const override
 Returns the AllocateContinuousState value, which must not be nullptr. More...
 
std::unique_ptr< DiscreteValues< T > > AllocateDiscreteVariables () const override
 Returns the AllocateDiscreteState value, which must not be nullptr. More...
 
std::multimap< int, intGetDirectFeedthroughs () const final
 Reports all direct feedthroughs from input ports to output ports. More...
 
int get_num_continuous_states () const final
 Returns the dimension of the continuous state vector that has been declared until now. More...
 
 LeafSystem (const LeafSystem &)=delete
 
LeafSystemoperator= (const LeafSystem &)=delete
 
 LeafSystem (LeafSystem &&)=delete
 
LeafSystemoperator= (LeafSystem &&)=delete
 
- Public Member Functions inherited from System< T >
 ~System () override=default
 
void GetWitnessFunctions (const Context< T > &context, std::vector< const WitnessFunction< T > *> *w) const
 Gets the witness functions active for the given state. More...
 
T CalcWitnessValue (const Context< T > &context, const WitnessFunction< T > &witness_func) const
 Evaluates a witness function at the given context. More...
 
 System (const System &)=delete
 
Systemoperator= (const System &)=delete
 
 System (System &&)=delete
 
Systemoperator= (System &&)=delete
 
std::unique_ptr< Context< T > > AllocateContext () const
 Returns a Context<T> suitable for use with this System<T>. More...
 
std::unique_ptr< BasicVector< T > > AllocateInputVector (const InputPort< T > &input_port) const
 Given an input port, allocates the vector storage. More...
 
std::unique_ptr< AbstractValueAllocateInputAbstract (const InputPort< T > &input_port) const
 Given an input port, allocates the abstract storage. More...
 
std::unique_ptr< SystemOutput< T > > AllocateOutput () const
 Returns a container that can hold the values of all of this System's output ports. More...
 
std::unique_ptr< Context< T > > CreateDefaultContext () const
 This convenience method allocates a context using AllocateContext() and sets its default values using SetDefaultContext(). More...
 
void SetDefaultContext (Context< T > *context) const
 
virtual void SetRandomState (const Context< T > &context, State< T > *state, RandomGenerator *generator) const
 Assigns random values to all elements of the state. More...
 
virtual void SetRandomParameters (const Context< T > &context, Parameters< T > *parameters, RandomGenerator *generator) const
 Assigns random values to all parameters. More...
 
void SetRandomContext (Context< T > *context, RandomGenerator *generator) const
 
void AllocateFixedInputs (Context< T > *context) const
 For each input port, allocates a fixed input of the concrete type that this System requires, and binds it to the port, disconnecting any prior input. More...
 
bool HasAnyDirectFeedthrough () const
 Returns true if any of the inputs to the system might be directly fed through to any of its outputs and false otherwise. More...
 
bool HasDirectFeedthrough (int output_port) const
 Returns true if there might be direct-feedthrough from any input port to the given output_port, and false otherwise. More...
 
bool HasDirectFeedthrough (int input_port, int output_port) const
 Returns true if there might be direct-feedthrough from the given input_port to the given output_port, and false otherwise. More...
 
void Publish (const Context< T > &context, const EventCollection< PublishEvent< T >> &events) const
 This method is the public entry point for dispatching all publish event handlers. More...
 
void Publish (const Context< T > &context) const
 Forces a publish on the system, given a context. More...
 
const ContinuousState< T > & EvalTimeDerivatives (const Context< T > &context) const
 Returns a reference to the cached value of the continuous state variable time derivatives, evaluating first if necessary using CalcTimeDerivatives(). More...
 
const TEvalPotentialEnergy (const Context< T > &context) const
 Returns a reference to the cached value of the potential energy (PE), evaluating first if necessary using CalcPotentialEnergy(). More...
 
const TEvalKineticEnergy (const Context< T > &context) const
 Returns a reference to the cached value of the kinetic energy (KE), evaluating first if necessary using CalcKineticEnergy(). More...
 
const TEvalConservativePower (const Context< T > &context) const
 Returns a reference to the cached value of the conservative power (Pc), evaluating first if necessary using CalcConservativePower(). More...
 
const TEvalNonConservativePower (const Context< T > &context) const
 Returns a reference to the cached value of the non-conservative power (Pnc), evaluating first if necessary using CalcNonConservativePower(). More...
 
template<template< typename > class Vec = BasicVector>
const Vec< T > * EvalVectorInput (const Context< T > &context, int port_index) const
 Returns the value of the vector-valued input port with the given port_index as a BasicVector or a specific subclass Vec derived from BasicVector. More...
 
Eigen::VectorBlock< const VectorX< T > > EvalEigenVectorInput (const Context< T > &context, int port_index) const
 Returns the value of the vector-valued input port with the given port_index as an Eigen vector. More...
 
int get_num_constraint_equations (const Context< T > &context) const
 Gets the number of constraint equations for this system using the given context (useful in case the number of constraints is dependent upon the current state (as might be the case with a system modeled using piecewise differential algebraic equations). More...
 
Eigen::VectorXd EvalConstraintEquations (const Context< T > &context) const
 Evaluates the constraint equations for the system at the generalized coordinates and generalized velocity specified by the context. More...
 
Eigen::VectorXd EvalConstraintEquationsDot (const Context< T > &context) const
 Computes the time derivative of each constraint equation, evaluated at the generalized coordinates and generalized velocity specified by the context. More...
 
Eigen::VectorXd CalcVelocityChangeFromConstraintImpulses (const Context< T > &context, const Eigen::MatrixXd &J, const Eigen::VectorXd &lambda) const
 Computes the change in velocity from applying the given constraint forces to the system at the given context. More...
 
double CalcConstraintErrorNorm (const Context< T > &context, const Eigen::VectorXd &error) const
 Computes the norm on constraint error (used as a metric for comparing errors between the outputs of algebraic equations applied to two different state variable instances). More...
 
void CalcTimeDerivatives (const Context< T > &context, ContinuousState< T > *derivatives) const
 Calculates the time derivatives xcdot of the continuous state xc into a given output argument. More...
 
void CalcDiscreteVariableUpdates (const Context< T > &context, const EventCollection< DiscreteUpdateEvent< T >> &events, DiscreteValues< T > *discrete_state) const
 This method is the public entry point for dispatching all discrete variable update event handlers. More...
 
void CalcDiscreteVariableUpdates (const Context< T > &context, DiscreteValues< T > *discrete_state) const
 This method forces a discrete update on the system given a context, and the updated discrete state is stored in discrete_state. More...
 
void CalcUnrestrictedUpdate (const Context< T > &context, const EventCollection< UnrestrictedUpdateEvent< T >> &events, State< T > *state) const
 This method is the public entry point for dispatching all unrestricted update event handlers. More...
 
void CalcUnrestrictedUpdate (const Context< T > &context, State< T > *state) const
 This method forces an unrestricted update on the system given a context, and the updated state is stored in state. More...
 
T CalcNextUpdateTime (const Context< T > &context, CompositeEventCollection< T > *events) const
 This method is called by a Simulator during its calculation of the size of the next continuous step to attempt. More...
 
void GetPerStepEvents (const Context< T > &context, CompositeEventCollection< T > *events) const
 This method is called by Simulator::Initialize() to gather all update and publish events that are to be handled in StepTo() at the point before Simulator integrates continuous state. More...
 
void GetInitializationEvents (const Context< T > &context, CompositeEventCollection< T > *events) const
 This method is called by Simulator::Initialize() to gather all update and publish events that need to be handled at initialization before the simulator starts integration. More...
 
optional< PeriodicEventDataGetUniquePeriodicDiscreteUpdateAttribute () const
 Gets whether there exists a unique periodic attribute that triggers one or more discrete update events (and, if so, returns that unique periodic attribute). More...
 
std::map< PeriodicEventData, std::vector< const Event< T > * >, PeriodicEventDataComparatorGetPeriodicEvents () const
 Gets all periodic triggered events for a system. More...
 
void CalcOutput (const Context< T > &context, SystemOutput< T > *outputs) const
 Utility method that computes for every output port i the value y(i) that should result from the current contents of the given Context. More...
 
T CalcPotentialEnergy (const Context< T > &context) const
 Calculates and returns the potential energy represented by the current configuration provided in context. More...
 
T CalcKineticEnergy (const Context< T > &context) const
 Calculates and returns the kinetic energy represented by the current configuration and velocity provided in context. More...
 
T CalcConservativePower (const Context< T > &context) const
 Calculates and returns the conservative power represented by the current contents of the given context. More...
 
T CalcNonConservativePower (const Context< T > &context) const
 Calculates and returns the non-conservative power represented by the current contents of the given context. More...
 
void MapVelocityToQDot (const Context< T > &context, const VectorBase< T > &generalized_velocity, VectorBase< T > *qdot) const
 Transforms a given generalized velocity v to the time derivative qdot of the generalized configuration q taken from the supplied Context. More...
 
void MapVelocityToQDot (const Context< T > &context, const Eigen::Ref< const VectorX< T >> &generalized_velocity, VectorBase< T > *qdot) const
 Transforms the given generalized velocity to the time derivative of generalized configuration. More...
 
void MapQDotToVelocity (const Context< T > &context, const VectorBase< T > &qdot, VectorBase< T > *generalized_velocity) const
 Transforms the time derivative qdot of the generalized configuration q to generalized velocities v. More...
 
void MapQDotToVelocity (const Context< T > &context, const Eigen::Ref< const VectorX< T >> &qdot, VectorBase< T > *generalized_velocity) const
 Transforms the given time derivative qdot of generalized configuration q to generalized velocity v. More...
 
std::string GetMemoryObjectName () const
 Returns a name for this System based on a stringification of its type name and memory address. More...
 
const InputPort< T > & get_input_port (int port_index) const
 Returns the typed input port at index port_index. More...
 
const InputPort< T > & GetInputPort (const std::string &port_name) const
 Returns the typed input port with the unique name port_name. More...
 
const OutputPort< T > & get_output_port (int port_index) const
 Returns the typed output port at index port_index. More...
 
const OutputPort< T > & GetOutputPort (const std::string &port_name) const
 Returns the typed output port with the unique name port_name. More...
 
int get_num_constraints () const
 Returns the number of constraints specified for the system. More...
 
const SystemConstraint< T > & get_constraint (SystemConstraintIndex constraint_index) const
 Returns the constraint at index constraint_index. More...
 
boolean< TCheckSystemConstraintsSatisfied (const Context< T > &context, double tol) const
 Returns true if context satisfies all of the registered SystemConstraints with tolerance tol. More...
 
void CheckValidOutput (const SystemOutput< T > *output) const
 Checks that output is consistent with the number and size of output ports declared by the system. More...
 
template<typename T1 = T>
void CheckValidContextT (const Context< T1 > &context) const
 Checks that context is consistent for this System template. More...
 
VectorX< TCopyContinuousStateVector (const Context< T > &context) const
 Returns a copy of the continuous state vector xc into an Eigen vector. More...
 
std::string GetGraphvizString (int max_depth=std::numeric_limits< int >::max()) const
 Returns a Graphviz string describing this System. More...
 
int64_t GetGraphvizId () const
 Returns an opaque integer that uniquely identifies this system in the Graphviz output. More...
 
void FixInputPortsFrom (const System< double > &other_system, const Context< double > &other_context, Context< T > *target_context) const
 Fixes all of the input ports in target_context to their current values in other_context, as evaluated by other_system. More...
 
const SystemScalarConverterget_system_scalar_converter () const
 (Advanced) Returns the SystemScalarConverter for this object. More...
 
std::unique_ptr< System< AutoDiffXd > > ToAutoDiffXd () const
 Creates a deep copy of this System, transmogrified to use the autodiff scalar type, with a dynamic-sized vector of partial derivatives. More...
 
std::unique_ptr< System< AutoDiffXd > > ToAutoDiffXdMaybe () const
 Creates a deep copy of this system exactly like ToAutoDiffXd(), but returns nullptr if this System does not support autodiff, instead of throwing an exception. More...
 
std::unique_ptr< System< symbolic::Expression > > ToSymbolic () const
 Creates a deep copy of this System, transmogrified to use the symbolic scalar type. More...
 
std::unique_ptr< System< symbolic::Expression > > ToSymbolicMaybe () const
 Creates a deep copy of this system exactly like ToSymbolic(), but returns nullptr if this System does not support symbolic, instead of throwing an exception. More...
 
- Public Member Functions inherited from SystemBase
 ~SystemBase () override
 
void set_name (const std::string &name)
 Sets the name of the system. More...
 
const std::string & get_name () const
 Returns the name last supplied to set_name(), if any. More...
 
const std::string & GetSystemName () const final
 Returns a human-readable name for this system, for use in messages and logging. More...
 
std::string GetSystemPathname () const final
 Generates and returns a human-readable full path name of this subsystem, for use in messages and logging. More...
 
std::string GetSystemType () const final
 Returns the most-derived type of this concrete System object as a human-readable string suitable for use in error messages. More...
 
void ThrowIfContextNotCompatible (const ContextBase &context) const final
 Throws an exception with an appropriate message if the given context is not compatible with this System. More...
 
std::unique_ptr< ContextBaseAllocateContext () const
 Returns a Context suitable for use with this System. More...
 
int get_num_input_ports () const
 Returns the number of input ports currently allocated in this System. More...
 
int get_num_output_ports () const
 Returns the number of output ports currently allocated in this System. More...
 
const InputPortBaseget_input_port_base (InputPortIndex port_index) const
 Returns a reference to an InputPort given its port_index. More...
 
const OutputPortBaseget_output_port_base (OutputPortIndex port_index) const
 Returns a reference to an OutputPort given its port_index. More...
 
int get_num_total_inputs () const
 Returns the total dimension of all of the vector-valued input ports (as if they were muxed). More...
 
int get_num_total_outputs () const
 Returns the total dimension of all of the vector-valued output ports (as if they were muxed). More...
 
int num_cache_entries () const
 Returns the number nc of cache entries currently allocated in this System. More...
 
const CacheEntryget_cache_entry (CacheIndex index) const
 Return a reference to a CacheEntry given its index. More...
 
void CheckValidContext (const ContextBase &context) const
 Checks whether the given context is valid for this System and throws an exception with a helpful message if not. More...
 
 SystemBase (const SystemBase &)=delete
 
SystemBaseoperator= (const SystemBase &)=delete
 
 SystemBase (SystemBase &&)=delete
 
SystemBaseoperator= (SystemBase &&)=delete
 
const AbstractValueEvalAbstractInput (const ContextBase &context, int port_index) const
 Returns the value of the input port with the given port_index as an AbstractValue, which is permitted for ports of any type. More...
 
template<typename V >
const V * EvalInputValue (const ContextBase &context, int port_index) const
 Returns the value of an abstract-valued input port with the given port_index as a value of known type V. More...
 
const CacheEntryDeclareCacheEntry (std::string description, CacheEntry::AllocCallback alloc_function, CacheEntry::CalcCallback calc_function, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a new CacheEntry in this System using the least-restrictive definitions for the associated functions. More...
 
template<class MySystem , class MyContext , typename ValueType >
const CacheEntryDeclareCacheEntry (std::string description, ValueType(MySystem::*make)() const, void(MySystem::*calc)(const MyContext &, ValueType *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a cache entry by specifying member functions to use both for the allocator and calculator. More...
 
template<class MySystem , class MyContext , typename ValueType >
const CacheEntryDeclareCacheEntry (std::string description, const ValueType &model_value, void(MySystem::*calc)(const MyContext &, ValueType *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a cache entry by specifying a model value of concrete type ValueType and a calculator function that is a class member function (method) with signature: More...
 
template<class MySystem , class MyContext , typename ValueType >
const CacheEntryDeclareCacheEntry (std::string description, const ValueType &model_value, ValueType(MySystem::*calc)(const MyContext &) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a cache entry by specifying a model value of concrete type ValueType and a calculator function that is a class member function (method) with signature: More...
 
template<class MySystem , class MyContext , typename ValueType >
const CacheEntryDeclareCacheEntry (std::string description, void(MySystem::*calc)(const MyContext &, ValueType *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a cache entry by specifying only a calculator function that is a class member function (method) with signature: More...
 
template<class MySystem , class MyContext , typename ValueType >
const CacheEntryDeclareCacheEntry (std::string description, ValueType(MySystem::*calc)(const MyContext &) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a cache entry by specifying only a calculator function that is a class member function (method) with signature: More...
 
DependencyTicket discrete_state_ticket (DiscreteStateIndex index) const
 Returns a ticket indicating dependence on a particular discrete state variable xdᵢ (may be a vector). More...
 
DependencyTicket abstract_state_ticket (AbstractStateIndex index) const
 Returns a ticket indicating dependence on a particular abstract state variable xaᵢ. More...
 
DependencyTicket numeric_parameter_ticket (NumericParameterIndex index) const
 Returns a ticket indicating dependence on a particular numeric parameter pnᵢ (may be a vector). More...
 
DependencyTicket abstract_parameter_ticket (AbstractParameterIndex index) const
 Returns a ticket indicating dependence on a particular abstract parameter paᵢ. More...
 
DependencyTicket input_port_ticket (InputPortIndex index)
 Returns a ticket indicating dependence on input port uᵢ indicated by index. More...
 
DependencyTicket cache_entry_ticket (CacheIndex index)
 Returns a ticket indicating dependence on the cache entry indicated by index. More...
 
DependencyTicket output_port_ticket (OutputPortIndex index)
 (Internal use only) Returns a ticket indicating dependence on the output port indicated by index. More...
 
int num_discrete_state_groups () const
 Returns the number of declared discrete state groups (each group is a vector-valued discrete state variable). More...
 
int num_abstract_states () const
 Returns the number of declared abstract state variables. More...
 
int num_numeric_parameters () const
 Returns the number of declared numeric parameters (each of these is a vector-valued parameter). More...
 
int num_abstract_parameters () const
 Returns the number of declared abstract parameters. More...
 

Friends

template<typename U >
class MultibodyPlant
 
class MultibodyPlantTester
 

Additional Inherited Members

- Static Public Member Functions inherited from System< T >
template<template< typename > class S = ::drake::systems::System>
static std::unique_ptr< S< AutoDiffXd > > ToAutoDiffXd (const S< T > &from)
 Creates a deep copy of from, transmogrified to use the autodiff scalar type, with a dynamic-sized vector of partial derivatives. More...
 
template<template< typename > class S = ::drake::systems::System>
static std::unique_ptr< S< symbolic::Expression > > ToSymbolic (const S< T > &from)
 Creates a deep copy of from, transmogrified to use the symbolic scalar type. More...
 
- Static Public Member Functions inherited from SystemBase
static DependencyTicket nothing_ticket ()
 Returns a ticket indicating that a computation does not depend on any source value; that is, it is a constant. More...
 
static DependencyTicket time_ticket ()
 Returns a ticket indicating dependence on time. More...
 
static DependencyTicket accuracy_ticket ()
 Returns a ticket indicating dependence on the accuracy setting in the Context. More...
 
static DependencyTicket q_ticket ()
 Returns a ticket indicating that a computation depends on configuration state variables q. More...
 
static DependencyTicket v_ticket ()
 Returns a ticket indicating dependence on velocity state variables v. More...
 
static DependencyTicket z_ticket ()
 Returns a ticket indicating dependence on any or all of the miscellaneous continuous state variables z. More...
 
static DependencyTicket xc_ticket ()
 Returns a ticket indicating dependence on all of the continuous state variables q, v, or z. More...
 
static DependencyTicket xd_ticket ()
 Returns a ticket indicating dependence on all of the numerical discrete state variables, in any discrete variable group. More...
 
static DependencyTicket xa_ticket ()
 Returns a ticket indicating dependence on all of the abstract state variables in the current Context. More...
 
static DependencyTicket all_state_ticket ()
 Returns a ticket indicating dependence on all state variables x in this system, including continuous variables xc, discrete (numeric) variables xd, and abstract state variables xa. More...
 
static DependencyTicket pn_ticket ()
 Returns a ticket indicating dependence on all of the numerical parameters in the current Context. More...
 
static DependencyTicket pa_ticket ()
 Returns a ticket indicating dependence on all of the abstract parameters pa in the current Context. More...
 
static DependencyTicket all_parameters_ticket ()
 Returns a ticket indicating dependence on all parameters p in this system, including numeric parameters pn, and abstract parameters pa. More...
 
static DependencyTicket all_input_ports_ticket ()
 Returns a ticket indicating dependence on all input ports u of this system. More...
 
static DependencyTicket all_sources_ticket ()
 Returns a ticket indicating dependence on every possible independent source value, including time, accuracy, state, input ports, and parameters (but not cache entries). More...
 
static DependencyTicket configuration_ticket ()
 Returns a ticket indicating dependence on all source values that may affect configuration-dependent computations. More...
 
static DependencyTicket kinematics_ticket ()
 Returns a ticket indicating dependence on all source values that may affect configuration- or velocity-dependent computations. More...
 
static DependencyTicket xcdot_ticket ()
 Returns a ticket for the cache entry that holds time derivatives of the continuous variables. More...
 
static DependencyTicket pe_ticket ()
 Returns a ticket for the cache entry that holds the potential energy calculation. More...
 
static DependencyTicket ke_ticket ()
 Returns a ticket for the cache entry that holds the kinetic energy calculation. More...
 
static DependencyTicket pc_ticket ()
 Returns a ticket for the cache entry that holds the conservative power calculation. More...
 
static DependencyTicket pnc_ticket ()
 Returns a ticket for the cache entry that holds the non-conservative power calculation. More...
 
- Protected Member Functions inherited from MultibodyTreeSystem< T >
void SetDefaultState (const systems::Context< T > &context, systems::State< T > *state) const override
 Assigns default values to all elements of the state. More...
 
 MultibodyTreeSystem (bool is_discrete=false)
 Default constructor allocates a MultibodyTree, with the intent that it will be filled in later, using mutable_tree() for access. More...
 
 MultibodyTreeSystem (systems::SystemScalarConverter converter, std::unique_ptr< MultibodyTree< T >> tree, bool is_discrete=false)
 Constructor that specifies scalar-type conversion support. More...
 
MultibodyTree< T > & mutable_tree () const
 Returns a mutable reference to the MultibodyTree owned by this class. More...
 
void Finalize ()
 Finalize the tree if that hasn't already been done, complete System construction, and declare any needed Context resources for the tree. More...
 
- Protected Member Functions inherited from LeafSystem< T >
 LeafSystem ()
 Default constructor that declares no inputs, outputs, state, parameters, events, nor scalar-type conversion support (AutoDiff, etc.). More...
 
 LeafSystem (SystemScalarConverter converter)
 Constructor that declares no inputs, outputs, state, parameters, or events, but allows subclasses to declare scalar-type conversion support (AutoDiff, etc.). More...
 
virtual void DoValidateAllocatedLeafContext (const LeafContext< T > &context) const
 Derived classes that impose restrictions on what resources are permitted should check those restrictions by implementing this. More...
 
T DoCalcWitnessValue (const Context< T > &context, const WitnessFunction< T > &witness_func) const final
 Derived classes will implement this method to evaluate a witness function at the given context. More...
 
void AddTriggeredWitnessFunctionToCompositeEventCollection (Event< T > *event, CompositeEventCollection< T > *events) const final
 Add event to events due to a witness function triggering. More...
 
void DoCalcNextUpdateTime (const Context< T > &context, CompositeEventCollection< T > *events, T *time) const override
 Computes the next update time based on the configured periodic events, for scalar types that are arithmetic, or aborts for scalar types that are not arithmetic. More...
 
BasicVector< T > * DoAllocateInputVector (const InputPort< T > &input_port) const override
 Allocates a vector that is suitable as an input value for input_port. More...
 
AbstractValueDoAllocateInputAbstract (const InputPort< T > &input_port) const override
 Allocates an AbstractValue suitable as an input value for input_port. More...
 
void GetGraphvizFragment (int max_depth, std::stringstream *dot) const override
 Emits a graphviz fragment for this System. More...
 
void GetGraphvizInputPortToken (const InputPort< T > &port, int max_depth, std::stringstream *dot) const final
 Appends a fragment to the dot stream identifying the graphviz node representing port. More...
 
void GetGraphvizOutputPortToken (const OutputPort< T > &port, int max_depth, std::stringstream *dot) const final
 Appends a fragment to the dot stream identifying the graphviz node representing port. More...
 
virtual std::unique_ptr< ContinuousState< T > > AllocateContinuousState () const
 Returns a ContinuousState used to implement both CreateDefaultContext and AllocateTimeDerivatives. More...
 
virtual std::unique_ptr< DiscreteValues< T > > AllocateDiscreteState () const
 Reserves the discrete state as required by CreateDefaultContext. More...
 
virtual std::unique_ptr< AbstractValuesAllocateAbstractState () const
 Reserves the abstract state as required by CreateDefaultContext. More...
 
virtual std::unique_ptr< Parameters< T > > AllocateParameters () const
 Reserves the parameters as required by CreateDefaultContext. More...
 
int DeclareNumericParameter (const BasicVector< T > &model_vector)
 Declares a numeric parameter using the given model_vector. More...
 
template<template< typename > class U = BasicVector>
const U< T > & GetNumericParameter (const Context< T > &context, int index) const
 Extracts the numeric parameters of type U from the context at index. More...
 
template<template< typename > class U = BasicVector>
U< T > & GetMutableNumericParameter (Context< T > *context, int index) const
 Extracts the numeric parameters of type U from the context at index. More...
 
int DeclareAbstractParameter (const AbstractValue &model_value)
 Declares an abstract parameter using the given model_value. More...
 
template<typename EventType >
void DeclarePeriodicEvent (double period_sec, double offset_sec)
 Declares that this System has a simple, fixed-period event specified with no custom callback function, and its attribute field contains an Event<T>::PeriodicAttribute constructed from the specified period_sec and offset_sec. More...
 
template<typename EventType >
void DeclarePeriodicEvent (double period_sec, double offset_sec, const EventType &event)
 Declares that this System has a simple, fixed-period event specified by event. More...
 
void DeclarePeriodicDiscreteUpdate (double period_sec, double offset_sec=0)
 Declares a periodic discrete update event with period = period_sec and offset = offset_sec. More...
 
void DeclarePeriodicUnrestrictedUpdate (double period_sec, double offset_sec=0)
 Declares a periodic unrestricted update event with period = period_sec and offset = offset_sec. More...
 
void DeclarePeriodicPublish (double period_sec, double offset_sec=0)
 Declares a periodic publish event with period = period_sec and offset = offset_sec. More...
 
template<typename EventType >
void DeclarePerStepEvent (const EventType &event)
 Declares a per-step event using event, which is deep copied (the copy is maintained by this). More...
 
template<typename EventType >
void DeclareInitializationEvent (const EventType &event)
 Declares an initialization event by deep copying event and storing it internally. More...
 
void DeclareContinuousState (int num_state_variables)
 Declares that this System should reserve continuous state with num_state_variables state variables, which have no second-order structure. More...
 
void DeclareContinuousState (int num_q, int num_v, int num_z)
 Declares that this System should reserve continuous state with num_q generalized positions, num_v generalized velocities, and num_z miscellaneous state variables. More...
 
void DeclareContinuousState (const BasicVector< T > &model_vector)
 Declares that this System should reserve continuous state with model_vector.size() miscellaneous state variables, stored in a vector Cloned from model_vector. More...
 
void DeclareContinuousState (const BasicVector< T > &model_vector, int num_q, int num_v, int num_z)
 Declares that this System should reserve continuous state with num_q generalized positions, num_v generalized velocities, and num_z miscellaneous state variables, stored in a vector Cloned from model_vector. More...
 
void DeclareContinuousState (std::unique_ptr< BasicVector< T >> model_vector, int num_q, int num_v, int num_z)
 Declares that this System should reserve continuous state with num_q generalized positions, num_v generalized velocities, and num_z miscellaneous state variables, stored in the a vector Cloned from model_vector. More...
 
void DeclareDiscreteState (int num_state_variables)
 Declares that this System should reserve discrete state with num_state_variables state variables. More...
 
int DeclareAbstractState (std::unique_ptr< AbstractValue > abstract_state)
 Declares an abstract state. More...
 
template<class MySystem >
SystemConstraintIndex DeclareEqualityConstraint (void(MySystem::*calc)(const Context< T > &, VectorX< T > *) const, int count, const std::string &description)
 Declares a system constraint of the form f(context) = 0 by specifying a member function to use to calculate the (VectorX) constraint value with a signature: More...
 
SystemConstraintIndex DeclareEqualityConstraint (typename SystemConstraint< T >::CalcCallback calc, int count, const std::string &description)
 Declares a system constraint of the form f(context) = 0 by specifying a std::function to use to calculate the (Vector) constraint value with a signature: More...
 
template<class MySystem >
SystemConstraintIndex DeclareInequalityConstraint (void(MySystem::*calc)(const Context< T > &, VectorX< T > *) const, int count, const std::string &description)
 Declares a system constraint of the form f(context) ≥ 0 by specifying a member function to use to calculate the (VectorX) constraint value with a signature: More...
 
SystemConstraintIndex DeclareInequalityConstraint (typename SystemConstraint< T >::CalcCallback calc, int count, const std::string &description)
 Declares a system constraint of the form f(context) ≥ 0 by specifying a std::function to use to calculate the (Vector) constraint value with a signature: More...
 
virtual void DoPublish (const Context< T > &context, const std::vector< const PublishEvent< T > *> &events) const
 Derived-class event handler for all simultaneous publish events in events. More...
 
virtual void DoCalcUnrestrictedUpdate (const Context< T > &context, const std::vector< const UnrestrictedUpdateEvent< T > *> &events, State< T > *state) const
 Derived-class event handler for all simultaneous unrestricted update events. More...
 
const InputPort< T > & DeclareVectorInputPort (std::string name, const BasicVector< T > &model_vector, optional< RandomDistribution > random_type=nullopt)
 Declares a vector-valued input port using the given model_vector. More...
 
const InputPort< T > & DeclareAbstractInputPort (std::string name, const AbstractValue &model_value)
 Declares an abstract-valued input port using the given model_value. More...
 
const InputPort< T > & DeclareVectorInputPort (const BasicVector< T > &model_vector, optional< RandomDistribution > random_type=nullopt)
 See the nearly identical signature with an additional (first) argument specifying the port name. More...
 
const InputPort< T > & DeclareAbstractInputPort (const AbstractValue &model_value)
 See the nearly identical signature with an additional (first) argument specifying the port name. More...
 
template<class MySystem , typename BasicVectorSubtype >
const OutputPort< T > & DeclareVectorOutputPort (std::string name, const BasicVectorSubtype &model_vector, void(MySystem::*calc)(const Context< T > &, BasicVectorSubtype *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a vector-valued output port by specifying (1) a model vector of type BasicVectorSubtype derived from BasicVector and initialized to the correct size and desired initial value, and (2) a calculator function that is a class member function (method) with signature: More...
 
template<class MySystem , typename BasicVectorSubtype >
const OutputPort< T > & DeclareVectorOutputPort (std::string name, void(MySystem::*calc)(const Context< T > &, BasicVectorSubtype *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares a vector-valued output port by specifying only a calculator function that is a class member function (method) with signature: More...
 
const OutputPort< T > & DeclareVectorOutputPort (std::string name, const BasicVector< T > &model_vector, typename LeafOutputPort< T >::CalcVectorCallback vector_calc_function, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 (Advanced) Declares a vector-valued output port using the given model_vector and a function for calculating the port's value at runtime. More...
 
template<class MySystem , typename OutputType >
const OutputPort< T > & DeclareAbstractOutputPort (std::string name, const OutputType &model_value, void(MySystem::*calc)(const Context< T > &, OutputType *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares an abstract-valued output port by specifying a model value of concrete type OutputType and a calculator function that is a class member function (method) with signature: More...
 
template<class MySystem , typename OutputType >
const OutputPort< T > & DeclareAbstractOutputPort (std::string name, void(MySystem::*calc)(const Context< T > &, OutputType *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares an abstract-valued output port by specifying only a calculator function that is a class member function (method) with signature: More...
 
template<class MySystem , typename OutputType >
const OutputPort< T > & DeclareAbstractOutputPort (std::string name, OutputType(MySystem::*make)() const, void(MySystem::*calc)(const Context< T > &, OutputType *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 Declares an abstract-valued output port by specifying member functions to use both for the allocator and calculator. More...
 
const OutputPort< T > & DeclareAbstractOutputPort (std::string name, typename LeafOutputPort< T >::AllocCallback alloc_function, typename LeafOutputPort< T >::CalcCallback calc_function, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 (Advanced) Declares an abstract-valued output port using the given allocator and calculator functions provided in their most generic forms. More...
 
template<class MySystem , typename BasicVectorSubtype >
const OutputPort< T > & DeclareVectorOutputPort (const BasicVectorSubtype &model_vector, void(MySystem::*calc)(const Context< T > &, BasicVectorSubtype *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 See the nearly identical signature with an additional (first) argument specifying the port name. More...
 
template<class MySystem , typename BasicVectorSubtype >
const OutputPort< T > & DeclareVectorOutputPort (void(MySystem::*calc)(const Context< T > &, BasicVectorSubtype *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 See the nearly identical signature with an additional (first) argument specifying the port name. More...
 
const OutputPort< T > & DeclareVectorOutputPort (const BasicVector< T > &model_vector, typename LeafOutputPort< T >::CalcVectorCallback vector_calc_function, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 See the nearly identical signature with an additional (first) argument specifying the port name. More...
 
template<class MySystem , typename OutputType >
std::enable_if_t<!std::is_same< OutputType, std::string >::value, const OutputPort< T > & > DeclareAbstractOutputPort (const OutputType &model_value, void(MySystem::*calc)(const Context< T > &, OutputType *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 See the nearly identical signature with an additional (first) argument specifying the port name. More...
 
template<class MySystem , typename OutputType >
const OutputPort< T > & DeclareAbstractOutputPort (void(MySystem::*calc)(const Context< T > &, OutputType *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 See the nearly identical signature with an additional (first) argument specifying the port name. More...
 
template<class MySystem , typename OutputType >
const OutputPort< T > & DeclareAbstractOutputPort (OutputType(MySystem::*make)() const, void(MySystem::*calc)(const Context< T > &, OutputType *) const, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 See the nearly identical signature with an additional (first) argument specifying the port name. More...
 
const OutputPort< T > & DeclareAbstractOutputPort (typename LeafOutputPort< T >::AllocCallback alloc_function, typename LeafOutputPort< T >::CalcCallback calc_function, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 See the nearly identical signature with an additional (first) argument specifying the port name. More...
 
template<class MySystem >
std::unique_ptr< WitnessFunction< T > > DeclareWitnessFunction (const std::string &description, const WitnessFunctionDirection &direction_type, T(MySystem::*calc)(const Context< T > &) const) const
 Constructs the witness function with the given description (used primarily for debugging and logging), direction type, and calculator function; and with no event object. More...
 
std::unique_ptr< WitnessFunction< T > > DeclareWitnessFunction (const std::string &description, const WitnessFunctionDirection &direction_type, std::function< T(const Context< T > &)> calc) const
 Constructs the witness function with the given description (used primarily for debugging and logging), direction type, and calculator function; and with no event object. More...
 
template<class MySystem >
std::unique_ptr< WitnessFunction< T > > DeclareWitnessFunction (const std::string &description, const WitnessFunctionDirection &direction_type, T(MySystem::*calc)(const Context< T > &) const, void(MySystem::*publish_callback)(const Context< T > &, const PublishEvent< T > &) const) const
 Constructs the witness function with the given description (used primarily for debugging and logging), direction type, calculator function, and publish event callback function for when this triggers. More...
 
template<class MySystem >
std::unique_ptr< WitnessFunction< T > > DeclareWitnessFunction (const std::string &description, const WitnessFunctionDirection &direction_type, T(MySystem::*calc)(const Context< T > &) const, void(MySystem::*du_callback)(const Context< T > &, const DiscreteUpdateEvent< T > &, DiscreteValues< T > *) const) const
 Constructs the witness function with the given description (used primarily for debugging and logging), direction type, calculator function, and discrete update event callback function for when this triggers. More...
 
template<class MySystem >
std::unique_ptr< WitnessFunction< T > > DeclareWitnessFunction (const std::string &description, const WitnessFunctionDirection &direction_type, T(MySystem::*calc)(const Context< T > &) const, void(MySystem::*uu_callback)(const Context< T > &, const UnrestrictedUpdateEvent< T > &, State< T > *) const) const
 Constructs the witness function with the given description (used primarily for debugging and logging), direction type, calculator function, and unrestricted update event callback function for when this triggers. More...
 
template<class MySystem >
std::unique_ptr< WitnessFunction< T > > DeclareWitnessFunction (const std::string &description, const WitnessFunctionDirection &direction_type, T(MySystem::*calc)(const Context< T > &) const, const Event< T > &e) const
 Constructs the witness function with the given description (used primarily for debugging and logging), direction type, and calculator function, and with an object corresponding to the event that is to be dispatched when this witness function triggers. More...
 
std::unique_ptr< WitnessFunction< T > > DeclareWitnessFunction (const std::string &description, const WitnessFunctionDirection &direction_type, std::function< T(const Context< T > &)> calc, const Event< T > &e) const
 Constructs the witness function with the given description (used primarily for debugging and logging), direction type, and calculator function, and with an object corresponding to the event that is to be dispatched when this witness function triggers. More...
 
- Protected Member Functions inherited from System< T >
virtual void DoGetWitnessFunctions (const Context< T > &, std::vector< const WitnessFunction< T > *> *) const
 Derived classes can override this method to provide witness functions active for the given state. More...
 
SystemConstraintIndex AddConstraint (std::unique_ptr< SystemConstraint< T >> constraint)
 Adds an already-created constraint to the list of constraints for this System. More...
 
const EventCollection< PublishEvent< T > > & get_forced_publish_events () const
 
const EventCollection< DiscreteUpdateEvent< T > > & get_forced_discrete_update_events () const
 
const EventCollection< UnrestrictedUpdateEvent< T > > & get_forced_unrestricted_update_events () const
 
void set_forced_publish_events (std::unique_ptr< EventCollection< PublishEvent< T >>> forced)
 
void set_forced_discrete_update_events (std::unique_ptr< EventCollection< DiscreteUpdateEvent< T >>> forced)
 
void set_forced_unrestricted_update_events (std::unique_ptr< EventCollection< UnrestrictedUpdateEvent< T >>> forced)
 
 System (SystemScalarConverter converter)
 Constructs an empty System base class object and allocates base class resources, possibly supporting scalar-type conversion support (AutoDiff, etc.) using converter. More...
 
const InputPort< T > & DeclareInputPort (std::string name, PortDataType type, int size, optional< RandomDistribution > random_type=nullopt)
 Adds a port with the specified type and size to the input topology. More...
 
const InputPort< T > & DeclareAbstractInputPort (std::string name)
 Adds an abstract-valued port to the input topology. More...
 
const InputPort< T > & DeclareInputPort (PortDataType type, int size, optional< RandomDistribution > random_type=nullopt)
 See the nearly identical signature with an additional (first) argument specifying the port name. More...
 
const InputPort< T > & DeclareAbstractInputPort ()
 See the nearly identical signature with an argument specifying the port name. More...
 
virtual T DoCalcPotentialEnergy (const Context< T > &context) const
 Override this method for physical systems to calculate the potential energy PE currently stored in the configuration provided in the given Context. More...
 
virtual T DoCalcKineticEnergy (const Context< T > &context) const
 Override this method for physical systems to calculate the kinetic energy KE currently present in the motion provided in the given Context. More...
 
virtual T DoCalcConservativePower (const Context< T > &context) const
 Override this method to return the rate Pc at which mechanical energy is being converted from potential energy to kinetic energy by this system in the given Context. More...
 
virtual T DoCalcNonConservativePower (const Context< T > &context) const
 Override this method to return the rate Pnc at which work W is done on the system by non-conservative forces. More...
 
virtual int do_get_num_constraint_equations (const Context< T > &context) const
 Gets the number of constraint equations for this system from the given context. More...
 
virtual Eigen::VectorXd DoEvalConstraintEquations (const Context< T > &context) const
 Evaluates the constraint equations for the system at the generalized coordinates and generalized velocity specified by the context. More...
 
virtual Eigen::VectorXd DoEvalConstraintEquationsDot (const Context< T > &context) const
 Computes the time derivative of each constraint equation, evaluated at the generalized coordinates and generalized velocity specified by the context. More...
 
virtual Eigen::VectorXd DoCalcVelocityChangeFromConstraintImpulses (const Context< T > &context, const Eigen::MatrixXd &J, const Eigen::VectorXd &lambda) const
 Computes the change in velocity from applying the given constraint forces to the system at the given context. More...
 
virtual double DoCalcConstraintErrorNorm (const Context< T > &context, const Eigen::VectorXd &error) const
 Computes the norm of the constraint error. More...
 
Eigen::VectorBlock< VectorX< T > > GetMutableOutputVector (SystemOutput< T > *output, int port_index) const
 Returns a mutable Eigen expression for a vector valued output port with index port_index in this system. More...
 
- Protected Member Functions inherited from SystemBase
 SystemBase ()=default
 (Internal use only) Default constructor. More...
 
void AddInputPort (std::unique_ptr< InputPortBase > port)
 (Internal use only) Adds an already-constructed input port to this System. More...
 
void AddOutputPort (std::unique_ptr< OutputPortBase > port)
 (Internal use only) Adds an already-constructed output port to this System. More...
 
std::string NextInputPortName (std::string given_name) const
 (Internal use only) Returns a name for the next input port, using the given name if it isn't kUseDefaultName, otherwise making up a name like "u3" from the next available input port index. More...
 
std::string NextOutputPortName (std::string given_name) const
 (Internal use only) Returns a name for the next output port, using the given name if it isn't kUseDefaultName, otherwise making up a name like "y3" from the next available output port index. More...
 
void AddDiscreteStateGroup (DiscreteStateIndex index)
 (Internal use only) Assigns a ticket to a new discrete variable group with the given index. More...
 
void AddAbstractState (AbstractStateIndex index)
 (Internal use only) Assigns a ticket to a new abstract state variable with the given index. More...
 
void AddNumericParameter (NumericParameterIndex index)
 (Internal use only) Assigns a ticket to a new numeric parameter with the given index. More...
 
void AddAbstractParameter (AbstractParameterIndex index)
 (Internal use only) Assigns a ticket to a new abstract parameter with the given index. More...
 
const CacheEntryDeclareCacheEntryWithKnownTicket (DependencyTicket known_ticket, std::string description, CacheEntry::AllocCallback alloc_function, CacheEntry::CalcCallback calc_function, std::set< DependencyTicket > prerequisites_of_calc={ all_sources_ticket()})
 (Internal use only) This is for cache entries associated with pre-defined tickets, for example the cache entry for time derivatives. More...
 
const internal::SystemParentServiceInterface * get_parent_service () const
 Returns a pointer to the service interface of the immediately enclosing Diagram if one has been set, otherwise nullptr. More...
 
DependencyTicket assign_next_dependency_ticket ()
 (Internal use only) Assigns the next unused dependency ticket number, unique only within a particular system. More...
 
const AbstractValueEvalAbstractInputImpl (const char *func, const ContextBase &context, InputPortIndex port_index) const
 (Internal use only) Shared code for updating an input port and returning a pointer to its abstract value, or nullptr if the port is not connected. More...
 
void ThrowNegativePortIndex (const char *func, int port_index) const
 Throws std::out_of_range to report a negative port_index that was passed to API method func. More...
 
void ThrowInputPortIndexOutOfRange (const char *func, InputPortIndex port_index) const
 Throws std::out_of_range to report bad input port_index that was passed to API method func. More...
 
void ThrowOutputPortIndexOutOfRange (const char *func, OutputPortIndex port_index) const
 Throws std::out_of_range to report bad output port_index that was passed to API method func. More...
 
void ThrowNotAVectorInputPort (const char *func, InputPortIndex port_index) const
 Throws std::logic_error because someone misused API method func, that is only allowed for declared-vector input ports, on an abstract port whose index is given here. More...
 
void ThrowInputPortHasWrongType (const char *func, InputPortIndex port_index, const std::string &expected_type, const std::string &actual_type) const
 Throws std::logic_error because someone called API method func claiming the input port had some value type that was wrong. More...
 
void ThrowCantEvaluateInputPort (const char *func, InputPortIndex port_index) const
 Throws std::logic_error because someone called API method func, that requires this input port to be evaluatable, but the port was neither fixed nor connected. More...
 
const InputPortBaseGetInputPortBaseOrThrow (const char *func, int port_index) const
 (Internal use only) Returns the InputPortBase at index port_index, throwing std::out_of_range we don't like the port index. More...
 
const OutputPortBaseGetOutputPortBaseOrThrow (const char *func, int port_index) const
 (Internal use only) Returns the OutputPortBase at index port_index, throwing std::out_of_range we don't like the port index. More...
 
void InitializeContextBase (ContextBase *context) const
 This method must be invoked from within derived class DoAllocateContext() implementations right after the concrete Context object has been allocated. More...
 
- Static Protected Member Functions inherited from SystemBase
static void set_parent_service (SystemBase *child, const internal::SystemParentServiceInterface *parent_service)
 (Internal use only) Declares that parent_service is the service interface of the Diagram that owns this subsystem. More...
 

Detailed Description

template<typename T>
class drake::multibody::multibody_plant::MultibodyPlant< T >

MultibodyPlant is a Drake system framework representation (see systems::System) for the model of a physical system consisting of a collection of interconnected bodies.

See Multibody Dynamics for an overview of concepts/notation.

MultibodyPlant provides a user-facing API to:

  • add bodies, joints, force elements, and constraints,
  • register geometries to a provided SceneGraph instance,
  • create and manipulate its Context,
  • perform Context-dependent computational queries.

System dynamics

The state of a multibody system x = [q; v] is given by its generalized positions vector q, of size nq (see num_positions()), and by its generalized velocities vector v, of size nv (see num_velocities()). As a Drake System, MultibodyPlant implements the governing equations for a multibody dynamical system in the form ẋ = f(t, x, u) with t being the time and u the input vector of actuation forces. The governing equations for the dynamics of a multibody system modeled with MultibodyPlant are [Featherstone 2008, Jain 2010]:

         q̇ = N(q)v
  (1)    M(q)v̇ + C(q, v)v = tau

where M(q) is the mass matrix of the multibody system, C(q, v)v corresponds to the bias term containing Coriolis and gyroscopic effects and N(q) is the kinematic coupling matrix describing the relationship between the rate of change of the generalized coordinates and the generalized velocities, [Seth 2010]. N(q) is an nq x nv matrix. The vector tau ∈ ℝⁿᵛ on the right hand side of Eq. (1) corresponds to generalized forces applied on the system. These can include externally applied body forces, constraint forces, and contact forces.

Loading models from SDF files

Drake has the capability of loading multibody models from SDF files. Consider the example below which loads an acrobot model from a file:

MultibodyPlant<T> acrobot;
const std::string relative_name =
"drake/multibody/benchmarks/acrobot/acrobot.sdf";
const std::string full_name = FindResourceOrThrow(relative_name);
AddModelFromSdfFile(full_name, &acrobot, &scene_graph);

As in the example above, for models including visual geometry, collision geometry or both, the user must specify a SceneGraph for geometry handling. You can find a full example of the LQR controlled acrobot in examples/multibody/acrobot/run_lqr.cc.

AddModelFromSdfFile() can be invoked multiple times on the same plant in order to load multiple model instances. Other parsing variants are available in multibody/multibody_tree/parsing/multibody_plant_sdf_parser.h such as AddModelsFromSdfFile() (please note the change to plural, i.e, "Models" instead of "Model") which allows creating model instances per each <model> tag found in the file. Please refer to each of these method's documentation for further details.

Adding modeling elements

Clients of a MultibodyPlant can add multibody elements with the following methods:

All modeling elements must be added pre-finalize.

Registering geometry with a SceneGraph

MultibodyPlant users can register geometry with a SceneGraph for essentially two purposes; a) visualization and, b) contact modeling. Before any geometry registration takes place, a user must first make a call to RegisterAsSourceForSceneGraph() in order to register the MultibodyPlant as a client of a SceneGraph instance, point at which the plant will have assigned a valid geometry::SourceId. At Finalize(), MultibodyPlant will declare input/output ports as appropriate to communicate with the SceneGraph instance on which registrations took place. All geometry registration must be performed pre-finalize.

If MultibodyPlant registers geometry with a SceneGraph via calls to RegisterCollisionGeometry(), an input port for geometric queries will be declared at Finalize() time, see get_geometry_query_input_port(). Users must connect this input port to the output port for geometric queries of the SceneGraph used for registration, which can be obtained with SceneGraph::get_query_output_port(). In summary, if MultibodyPlant registers collision geometry, the setup process will include:

  1. Call to RegisterAsSourceForSceneGraph().
  2. Calls to RegisterCollisionGeometry(), as many as needed.
  3. Call to Finalize(), user is done specifying the model.
  4. Connect SceneGraph::get_query_output_port() to get_geometry_query_input_port(). Refer to the documentation provided in each of the methods above for further details.

Constructor & Destructor Documentation

◆ MultibodyPlant() [1/4]

MultibodyPlant ( const MultibodyPlant< T > &  )
delete

◆ MultibodyPlant() [2/4]

MultibodyPlant ( MultibodyPlant< T > &&  )
delete

◆ MultibodyPlant() [3/4]

MultibodyPlant ( double  time_step = 0)
explicit

Default constructor creates a plant with a single "world" body.

Therefore, right after creation, num_bodies() returns one.

Parameters
[in]time_stepAn optional parameter indicating whether this plant is modeled as a continuous system (time_step = 0) or as a discrete system with periodic updates of period time_step > 0.
Default: 0.0.
Exceptions
std::exceptionif time_step is negative.

◆ MultibodyPlant() [4/4]

MultibodyPlant ( const MultibodyPlant< U > &  other)
inline

Scalar-converting copy constructor. See System Scalar Conversion.

Member Function Documentation

◆ AddForceElement() [1/2]

const ForceElementType<T>& AddForceElement ( Args &&...  args)
inline

Adds a new force element model of type ForceElementType to this MultibodyPlant.

The arguments to this method args are forwarded to ForceElementType's constructor.

Parameters
[in]argsZero or more parameters provided to the constructor of the new force element. It must be the case that JointType<T>(args) is a valid constructor.
Template Parameters
ForceElementTypeThe type of the ForceElement to add. This method can only be called once for elements of type UniformGravityFieldElement. That is, gravity can only be specified once.
Returns
A constant reference to the new ForceElement just added, of type ForceElementType<T> specialized on the scalar type T of this MultibodyPlant. It will remain valid for the lifetime of this MultibodyPlant.
See also
The ForceElement class's documentation for further details on how a force element is defined.

◆ AddForceElement() [2/2]

std::enable_if< std::is_same<ForceElementType<T>, UniformGravityFieldElement<T> >::value, const ForceElementType<T>&>::type AddForceElement ( Args &&...  args)
inline

◆ AddFrame()

const FrameType<T>& AddFrame ( std::unique_ptr< FrameType< T >>  frame)
inline

This method adds a Frame of type FrameType<T>.

For more information, please see the corresponding constructor of FrameType.

Template Parameters
FrameTypeTemplate which will be instantiated on T.
Parameters
frameUnique pointer frame instance.
Returns
A constant reference to the new Frame just added, which will remain valid for the lifetime of this MultibodyPlant.

◆ AddJoint() [1/2]

const JointType<T>& AddJoint ( std::unique_ptr< JointType< T >>  joint)
inline

This method adds a Joint of type JointType between two bodies.

For more information, see the below overload of AddJoint<>, and the related MultibodyTree::AddJoint<> method.

◆ AddJoint() [2/2]

const JointType<T>& AddJoint ( const std::string &  name,
const Body< T > &  parent,
const optional< Isometry3< double >> &  X_PF,
const Body< T > &  child,
const optional< Isometry3< double >> &  X_BM,
Args &&...  args 
)
inline

This method adds a Joint of type JointType between two bodies.

The two bodies connected by this Joint object are referred to as the parent and child bodies. Although the terms parent and child are sometimes used synonymously to describe the relationship between inboard and outboard bodies in multibody models, this usage is wholly unrelated and implies nothing about the inboard-outboard relationship between the bodies. As explained in the Joint class's documentation, in Drake we define a frame F attached to the parent body P with pose X_PF and a frame M attached to the child body B with pose X_BM. This method helps creating a joint between two bodies with fixed poses X_PF and X_BM. Refer to the Joint class's documentation for more details.

Parameters
nameA string that uniquely identifies the new joint to be added to this model. A std::runtime_error is thrown if a joint named name already is part of the model. See HasJointNamed(), Joint::name().
[in]parentThe parent body connected by the new joint.
[in]X_PFThe fixed pose of frame F attached to the parent body, measured in the frame P of that body. X_PF is an optional parameter; empty curly braces {} imply that frame F is the same body frame P. If instead your intention is to make a frame F with pose X_PF equal to the identity pose, provide Isometry3<double>::Identity() as your input.
[in]childThe child body connected by the new joint.
[in]X_BMThe fixed pose of frame M attached to the child body, measured in the frame B of that body. X_BM is an optional parameter; empty curly braces {} imply that frame M is the same body frame B. If instead your intention is to make a frame M with pose X_BM equal to the identity pose, provide Isometry3<double>::Identity() as your input.
[in]argsZero or more parameters provided to the constructor of the new joint. It must be the case that JointType<T>( const std::string&, const Frame<T>&, const Frame<T>&, args) is a valid constructor.
Template Parameters
JointTypeThe type of the Joint to add.
Returns
A constant reference to the new joint just added, of type JointType<T> specialized on the scalar type T of this MultibodyPlant. It will remain valid for the lifetime of this MultibodyPlant.

Example of usage:

MultibodyPlant<T> plant;
// Code to define bodies serving as the joint's parent and child bodies.
const RigidBody<double>& body_1 =
plant.AddRigidBody("Body1", SpatialInertia<double>(...));
const RigidBody<double>& body_2 =
plant.AddRigidBody("Body2", SpatialInertia<double>(...));
// Body 1 serves as parent, Body 2 serves as child.
// Define the pose X_BM of a frame M rigidly atached to child body B.
const RevoluteJoint<double>& elbow =
plant.AddJoint<RevoluteJoint>(
"Elbow", /* joint name */
body_1, /* parent body */
{}, /* frame F IS the parent body frame P */
body_2, /* child body, the pendulum */
X_BM, /* pose of frame M in the body frame B */
Vector3d::UnitZ()); /* revolute axis in this case */
Exceptions
std::exceptionif this MultibodyPlant already contains a joint with the given name. See HasJointNamed(), Joint::name().
See also
The Joint class's documentation for further details on how a Joint is defined.

◆ AddJointActuator()

const JointActuator<T>& AddJointActuator ( const std::string &  name,
const Joint< T > &  joint 
)
inline

Creates and adds a JointActuator model for an actuator acting on a given joint.

This method returns a constant reference to the actuator just added, which will remain valid for the lifetime of this plant.

Parameters
[in]nameA string that uniquely identifies the new actuator to be added to this model. A std::runtime_error is thrown if an actuator with the same name already exists in the model. See HasJointActuatorNamed().
[in]jointThe Joint to be actuated by the new JointActuator.
Returns
A constant reference to the new JointActuator just added, which will remain valid for the lifetime of this plant.
Exceptions
std::exceptionif joint.num_velocities() > 1 since for now we only support actuators for single dof joints.

◆ AddModelInstance()

ModelInstanceIndex AddModelInstance ( const std::string &  name)
inline

Creates a new model instance.

Returns the index for the model instance.

Parameters
[in]nameA string that uniquely identifies the new instance to be added to this model. An exception is thrown if an instance with the same name already exists in the model. See HasModelInstanceNamed().

◆ AddRigidBody() [1/2]

const RigidBody<T>& AddRigidBody ( const std::string &  name,
ModelInstanceIndex  model_instance,
const SpatialInertia< double > &  M_BBo_B 
)
inline

Creates a rigid body with the provided name and spatial inertia.

This method returns a constant reference to the body just added, which will remain valid for the lifetime of this MultibodyPlant.

Example of usage:

MultibodyPlant<T> plant;
// ... Code to define spatial_inertia, a SpatialInertia<T> object ...
ModelInstanceIndex model_instance = plant.AddModelInstance("instance");
const RigidBody<T>& body =
plant.AddRigidBody("BodyName", model_instance, spatial_inertia);
Parameters
[in]nameA string that identifies the new body to be added to this model. A std::runtime_error is thrown if a body named name already is part of model_instance. See HasBodyNamed(), Body::name().
[in]model_instanceA model instance index which this body is part of.
[in]M_BBo_BThe SpatialInertia of the new rigid body to be added to this MultibodyPlant, computed about the body frame origin Bo and expressed in the body frame B.
Returns
A constant reference to the new RigidBody just added, which will remain valid for the lifetime of this MultibodyPlant.

◆ AddRigidBody() [2/2]

const RigidBody<T>& AddRigidBody ( const std::string &  name,
const SpatialInertia< double > &  M_BBo_B 
)
inline

Creates a rigid body with the provided name and spatial inertia.

This method returns a constant reference to the body just added, which will remain valid for the lifetime of this MultibodyPlant. The body will use the default model instance (more on model instances).

Example of usage:

MultibodyPlant<T> plant;
// ... Code to define spatial_inertia, a SpatialInertia<T> object ...
const RigidBody<T>& body =
plant.AddRigidBody("BodyName", spatial_inertia);
Parameters
[in]nameA string that identifies the new body to be added to this model. A std::runtime_error is thrown if a body named name already is part of the model in the default model instance. See HasBodyNamed(), Body::name().
[in]M_BBo_BThe SpatialInertia of the new rigid body to be added to this MultibodyPlant, computed about the body frame origin Bo and expressed in the body frame B.
Returns
A constant reference to the new RigidBody just added, which will remain valid for the lifetime of this MultibodyPlant.
Exceptions
std::logic_errorif additional model instances have been created beyond the world and default instances.

◆ CollectRegisteredGeometries()

geometry::GeometrySet CollectRegisteredGeometries ( const std::vector< const Body< T > *> &  bodies) const

For each of the provided bodies, collects up all geometries that have been registered to that body.

Intended to be used in conjunction with SceneGraph::ExcludeCollisionsWithin() and SceneGraph::ExcludeCollisionsBetween() to filter collisions between the geometries registered to the bodies.

For example:

// Don't report on collisions between geometries affixed to `body1`,
// `body2`, or `body3`.
std::vector<const RigidBody<T>*> bodies{&body1, &body2, &body3};
geometry::GeometrySet set = plant.CollectRegisteredGeometries(bodies);
scene_graph.ExcludeCollisionsWithin(set);

Note: There is a very specific order of operations.

  1. Bodies and geometries must be added to the MultibodyPlant.
  2. The MultibodyPlant must be finalized (via Finalize()).
  3. Create GeometrySet instances from bodies (via this method).
  4. Invoke SceneGraph::ExcludeCollisions*() to filter collisions.
  5. Allocate context. Changing the order will cause exceptions to be thrown.
Exceptions
std::exceptionif called pre-finalize.

◆ default_coulomb_friction()

const CoulombFriction<double>& default_coulomb_friction ( geometry::GeometryId  id) const
inline

Returns the friction coefficients provided during geometry registration for the given geometry id.

We call these the "default" coefficients but note that we mean user-supplied per-geometry default, not something more global.

Exceptions
std::exceptionif id does not correspond to a geometry in this model registered for contact modeling.
See also
RegisterCollisionGeometry() for details on geometry registration.

◆ Finalize()

void Finalize ( geometry::SceneGraph< T > *  scene_graph = nullptr)

This method must be called after all elements in the model (joints, bodies, force elements, constraints, etc.) are added and before any computations are performed.

It essentially compiles all the necessary "topological information", i.e. how bodies, joints and, any other elements connect with each other, and performs all the required pre-processing to enable computations at a later stage.

If the finalize stage is successful, the topology of this MultibodyPlant is valid, meaning that the topology is up-to-date after this call. No more multibody elements can be added after a call to Finalize().

At Finalize(), state and input/output ports for this plant are declared. If this plant registered geometry with a SceneGraph, input and output ports to enable communication with that SceneGraph are declared as well.

If geometry has been registered on a SceneGraph instance, that instance must be provided to the Finalize() method so that any geometric implications of the finalization process can be appropriately handled.

See also
is_finalized().
Exceptions
std::logic_errorif
  1. the MultibodyPlant has already been finalized,
  2. scene_graph isn't provided when required, or
  3. a different scene_graph instance is provided than the one for which this plant is a geometry source.

◆ geometry_source_is_registered()

bool geometry_source_is_registered ( ) const
inline

Returns true if this MultibodyPlant was registered with a SceneGraph.

This method can be called at any time during the lifetime of this plant to query if this plant has been registered with a SceneGraph, either pre- or post-finalize, see Finalize().

◆ get_actuation_input_port() [1/2]

const systems::InputPort< T > & get_actuation_input_port ( ) const

Returns a constant reference to the input port for external actuation for the case where only one model instance has actuated dofs.

This input port is a vector valued port, which can be set with JointActuator::set_actuation_vector().

Precondition
Finalize() was already called on this plant.
Exceptions
std::exceptionif called before Finalize(), if the model does not contain any actuators, or if multiple model instances have actuated dofs.

◆ get_actuation_input_port() [2/2]

const systems::InputPort< T > & get_actuation_input_port ( ModelInstanceIndex  model_instance) const

Returns a constant reference to the input port for external actuation for a specific model instance.

This input port is a vector valued port, which can be set with JointActuator::set_actuation_vector().

Precondition
Finalize() was already called on this plant.
Exceptions
std::exceptionif called before Finalize() or if the model instance does not contain any actuators.
std::exceptionif the model instance does not exist.

◆ get_contact_penalty_method_time_scale()

double get_contact_penalty_method_time_scale ( ) const
inline

Returns a time-scale estimate tc based on the requested penetration allowance δ set with set_penetration_allowance().

For the penalty method in use to enforce non-penetration, this time scale relates to the time it takes the relative normal velocity between two bodies to go to zero. This time scale tc is artificially introduced by the penalty method and goes to zero in the limit to ideal rigid contact. Since numerical integration methods for continuum systems must be able to resolve a system's dynamics, the time step used by an integrator must in general be much smaller than the time scale tc. How much smaller will depend on the details of the problem and the convergence characteristics of the integrator and should be tuned appropriately. Another factor to take into account for setting up the simulation's time step is the speed of the objects in your simulation. If vn represents a reference velocity scale for the normal relative velocity between bodies, the new time scale tn = δ / vn represents the time it would take for the distance between two bodies approaching with relative normal velocity vn to decrease by the penetration_allowance δ. In this case a user should choose a time step for simulation that can resolve the smallest of the two time scales tc and tn.

◆ get_contact_results_output_port()

const systems::OutputPort< T > & get_contact_results_output_port ( ) const

Returns a constant reference to the port that outputs ContactResults.

Exceptions
std::exceptionif this plant is not modeled as a discrete system with periodic updates.
std::exceptionif called pre-finalize, see Finalize().

◆ get_continuous_state_output_port() [1/2]

const systems::OutputPort< T > & get_continuous_state_output_port ( ) const

Returns a constant reference to the output port for the full continuous state of the model.

Precondition
Finalize() was already called on this plant.

◆ get_continuous_state_output_port() [2/2]

const systems::OutputPort< T > & get_continuous_state_output_port ( ModelInstanceIndex  model_instance) const

Returns a constant reference to the output port for the continuous state of a specific model instance.

Precondition
Finalize() was already called on this plant.
Exceptions
std::exceptionif called before Finalize() or if the model instance does not have any state.
std::exceptionif the model instance does not exist.

◆ get_generalized_contact_forces_output_port()

const systems::OutputPort< T > & get_generalized_contact_forces_output_port ( ModelInstanceIndex  model_instance) const

Returns a constant reference to the output port of generalized contact forces for a specific model instance.

This output port is only available when modeling the plant as a discrete system with periodic updates, see is_discrete().

Precondition
Finalize() was already called on this plant.
Exceptions
std::exceptionif this plant is not modeled as a discrete system with periodic updates.
std::exceptionif called before Finalize() or if the model instance does not have any generalized velocities.
std::exceptionif the model instance does not exist.

◆ get_geometry_poses_output_port()

const OutputPort< T > & get_geometry_poses_output_port ( ) const

Returns the output port of frames' poses to communicate with a SceneGraph.

Exceptions
std::exceptionif this system was not registered with a SceneGraph.
std::exceptionif called pre-finalize. See Finalize().

◆ get_geometry_query_input_port()

const systems::InputPort< T > & get_geometry_query_input_port ( ) const

Returns a constant reference to the input port used to perform geometric queries on a SceneGraph.

See SceneGraph::get_query_output_port(). Refer to section Registering geometry with a SceneGraph of this class's documentation for further details on collision geometry registration and connection with a SceneGraph.

Exceptions
std::exceptionif this system was not registered with a SceneGraph.
std::exceptionif called pre-finalize. See Finalize().

◆ get_source_id()

optional<geometry::SourceId> get_source_id ( ) const
inline

Returns the unique id identifying this plant as a source for a SceneGraph.

Returns nullopt if this plant did not register any geometry. This method can be called at any time during the lifetime of this plant to query if this plant has been registered with a SceneGraph, either pre- or post-finalize, see Finalize(). However, a geometry::SourceId is only assigned once at the first call of any of this plant's geometry registration methods, and it does not change after that. Post-finalize calls will always return the same value.

◆ GetBodyByName() [1/2]

const Body<T>& GetBodyByName ( const std::string &  name) const
inline

These queries can be performed at any time during the lifetime of a MultibodyPlant, i.e.

there is no restriction on whether they must be called before or after Finalize(). This implies that these queries can be performed while new multibody elements are being added to the model.

If the named element is present in more than one model instance and a model instance is not explicitly specified, std::logic_error is thrown.Returns a constant reference to a body that is identified by the string name in this MultibodyPlant.

Exceptions
std::logic_errorif there is no body with the requested name.
std::logic_errorif the body name occurs in multiple model instances.
See also
HasBodyNamed() to query if there exists a body in this MultibodyPlant with a given specified name.

◆ GetBodyByName() [2/2]

const Body<T>& GetBodyByName ( const std::string &  name,
ModelInstanceIndex  model_instance 
) const
inline

Returns a constant reference to the body that is uniquely identified by the string name and model_instance in this MultibodyPlant.

Exceptions
std::logic_errorif there is no body with the requested name.
See also
HasBodyNamed() to query if there exists a body in this MultibodyPlant with a given specified name.

◆ GetBodyFrameIdIfExists()

optional<geometry::FrameId> GetBodyFrameIdIfExists ( BodyIndex  body_index) const
inline

If the body with body_index has geometry registered with it, it returns the geometry::FrameId associated with it.

Otherwise, it returns nullopt.

Exceptions
std::exceptionif called pre-finalize.

◆ GetBodyFrameIdOrThrow()

geometry::FrameId GetBodyFrameIdOrThrow ( BodyIndex  body_index) const
inline

If the body with body_index has geometry registered with it, it returns the geometry::FrameId associated with it.

Otherwise this method throws an exception.

Exceptions
std::exceptionif no geometry has been registered with the body indicated by body_index.
std::exceptionif called pre-finalize.

◆ GetBodyFromFrameId()

const Body<T>* GetBodyFromFrameId ( geometry::FrameId  frame_id) const
inline

Given a geometry frame identifier, returns a pointer to the body associated with that id (nullptr if there is no such body).

◆ GetCollisionGeometriesForBody()

const std::vector< geometry::GeometryId > & GetCollisionGeometriesForBody ( const Body< T > &  body) const

Returns an array of GeometryId's identifying the different contact geometries for body previously registered with a SceneGraph.

Note
This method can be called at any time during the lifetime of this plant, either pre- or post-finalize, see Finalize(). Post-finalize calls will always return the same value.
See also
RegisterCollisionGeometry(), Finalize()

◆ GetFrameByName() [1/2]

const Frame<T>& GetFrameByName ( const std::string &  name) const
inline

Returns a constant reference to a frame that is identified by the string name in this model.

Exceptions
std::logic_errorif there is no frame with the requested name.
std::logic_errorif the frame name occurs in multiple model instances.
See also
HasFrameNamed() to query if there exists a frame in this model with a given specified name.

◆ GetFrameByName() [2/2]

const Frame<T>& GetFrameByName ( const std::string &  name,
ModelInstanceIndex  model_instance 
) const
inline

Returns a constant reference to the frame that is uniquely identified by the string name in model_instance.

Exceptions
std::logic_errorif there is no frame with the requested name.
std::runtime_errorif model_instance is not valid for this model.
See also
HasFrameNamed() to query if there exists a frame in this model with a given specified name.

◆ GetJointActuatorByName() [1/2]

const JointActuator<T>& GetJointActuatorByName ( const std::string &  name) const
inline

Returns a constant reference to an actuator that is identified by the string name in this MultibodyPlant.

Exceptions
std::logic_errorif there is no actuator with the requested name.
std::logic_errorif the actuator name occurs in multiple model instances.
See also
HasJointActuatorNamed() to query if there exists an actuator in this MultibodyPlant with a given specified name.

◆ GetJointActuatorByName() [2/2]

const JointActuator<T>& GetJointActuatorByName ( const std::string &  name,
ModelInstanceIndex  model_instance 
) const
inline

Returns a constant reference to the actuator that is uniquely identified by the string name and model_instance in this MultibodyPlant.

Exceptions
std::logic_errorif there is no actuator with the requested name.
std::exceptionif model_instance is not valid for this model.
See also
HasJointActuatorNamed() to query if there exists an actuator in this MultibodyPlant with a given specified name.

◆ GetJointByName() [1/4]

const Joint<T>& GetJointByName ( const std::string &  name) const
inline

Returns a constant reference to a joint that is identified by the string name in this MultibodyPlant.

Exceptions
std::logic_errorif there is no joint with the requested name.
std::logic_errorif the joint name occurs in multiple model instances.
See also
HasJointNamed() to query if there exists a joint in this MultibodyPlant with a given specified name.

◆ GetJointByName() [2/4]

const Joint<T>& GetJointByName ( const std::string &  name,
ModelInstanceIndex  model_instance 
) const
inline

Returns a constant reference to the joint that is uniquely identified by the string name and model_instance in this MultibodyPlant.

Exceptions
std::logic_errorif there is no joint with the requested name.
std::exceptionif model_instance is not valid for this model.
See also
HasJointNamed() to query if there exists a joint in this MultibodyPlant with a given specified name.

◆ GetJointByName() [3/4]

const JointType<T>& GetJointByName ( const std::string &  name) const
inline

A templated version of GetJointByName() to return a constant reference of the specified type JointType in place of the base Joint class.

See GetJointByName() for details.

Template Parameters
JointTypeThe specific type of the Joint to be retrieved. It must be a subclass of Joint.
Exceptions
std::logic_errorif the named joint is not of type JointType or if there is no Joint with that name.
std::logic_errorif the joint name occurs in multiple model instances.
See also
HasJointNamed() to query if there exists a joint in this MultibodyPlant with a given specified name.

◆ GetJointByName() [4/4]

const JointType<T>& GetJointByName ( const std::string &  name,
ModelInstanceIndex  model_instance 
) const
inline

A templated version of GetJointByName() to return a constant reference of the specified type JointType in place of the base Joint class.

See GetJointByName() for details.

Template Parameters
JointTypeThe specific type of the Joint to be retrieved. It must be a subclass of Joint.
Exceptions
std::logic_errorif the named joint is not of type JointType or if there is no Joint with that name.
std::exceptionif model_instance is not valid for this model.
See also
HasJointNamed() to query if there exists a joint in this MultibodyPlant with a given specified name.

◆ GetModelInstanceByName()

ModelInstanceIndex GetModelInstanceByName ( const std::string &  name) const
inline

Returns the index to the model instance that is uniquely identified by the string name in this MultibodyPlant.

Exceptions
std::logic_errorif there is no instance with the requested name.
See also
HasModelInstanceNamed() to query if there exists an instance in this MultibodyPlant with a given specified name.

◆ GetVisualGeometriesForBody()

const std::vector< geometry::GeometryId > & GetVisualGeometriesForBody ( const Body< T > &  body) const

Returns an array of GeometryId's identifying the different visual geometries for body previously registered with a SceneGraph.

Note
This method can be called at any time during the lifetime of this plant, either pre- or post-finalize, see Finalize(). Post-finalize calls will always return the same value.
See also
RegisterVisualGeometry(), Finalize()

◆ HasBodyNamed() [1/2]

bool HasBodyNamed ( const std::string &  name) const
inline
Returns
true if a body named name was added to the MultibodyPlant.
See also
AddRigidBody().
Exceptions
std::logic_errorif the body name occurs in multiple model instances.

◆ HasBodyNamed() [2/2]

bool HasBodyNamed ( const std::string &  name,
ModelInstanceIndex  model_instance 
) const
inline
Returns
true if a body named name was added to the MultibodyPlant in model_instance.
See also
AddRigidBody().
Exceptions
std::exceptionif model_instance is not valid for this model.

◆ HasJointActuatorNamed() [1/2]

bool HasJointActuatorNamed ( const std::string &  name) const
inline
Returns
true if an actuator named name was added to the MultibodyPlant.
See also
AddJointActuator().
Exceptions
std::logic_errorif the actuator name occurs in multiple model instances.

◆ HasJointActuatorNamed() [2/2]

bool HasJointActuatorNamed ( const std::string &  name,
ModelInstanceIndex  model_instance 
) const
inline
Returns
true if an actuator named name was added to the MultibodyPlant in model_instance.
See also
AddJointActuator().
Exceptions
std::exceptionif model_instance is not valid for this model.

◆ HasJointNamed() [1/2]

bool HasJointNamed ( const std::string &  name) const
inline
Returns
true if a joint named name was added to the MultibodyPlant.
See also
AddJoint().
Exceptions
std::logic_errorif the joint name occurs in multiple model instances.

◆ HasJointNamed() [2/2]

bool HasJointNamed ( const std::string &  name,
ModelInstanceIndex  model_instance 
) const
inline
Returns
true if a joint named name was added to the MultibodyPlant in model_instance.
See also
AddJoint().
Exceptions
std::exceptionif model_instance is not valid for this model.

◆ HasModelInstanceNamed()

bool HasModelInstanceNamed ( const std::string &  name) const
inline
Returns
true if a model instance named name was added to the MultibodyPlant.
See also
AddModelInstance().

◆ is_finalized()

bool is_finalized ( ) const
inline

Returns true if this MultibodyPlant was finalized with a call to Finalize().

See also
Finalize().

◆ model()

const MultibodyTree<T>& model ( ) const
inline

Returns a constant reference to the underlying MultibodyTree model for this plant.

Exceptions
std::exceptionif called pre-finalize. See Finalize().

◆ num_actuated_dofs() [1/2]

int num_actuated_dofs ( ) const
inline

Returns the total number of actuated degrees of freedom.

That is, the vector of actuation values u has this size. See AddJointActuator().

◆ num_actuated_dofs() [2/2]

int num_actuated_dofs ( ModelInstanceIndex  model_instance) const
inline

Returns the total number of actuated degrees of freedom for a specific model instance.

That is, the vector of actuation values u has this size. See AddJointActuator().

◆ num_actuators()

int num_actuators ( ) const
inline

Returns the number of joint actuators in the model.

See also
AddJointActuator().

◆ num_bodies()

int num_bodies ( ) const
inline

Returns the number of bodies in the model, including the "world" body, which is always part of the model.

See also
AddRigidBody().

◆ num_collision_geometries()

int num_collision_geometries ( ) const
inline

Returns the number of geometries registered for contact modeling.

This method can be called at any time during the lifetime of this plant, either pre- or post-finalize, see Finalize(). Post-finalize calls will always return the same value.

◆ num_joints()

int num_joints ( ) const
inline

Returns the number of joints in the model.

See also
AddJoint().

◆ num_model_instances()

int num_model_instances ( ) const
inline

Returns the number of model instances in the model.

See also
AddModelInstance().

◆ num_multibody_states()

int num_multibody_states ( ) const
inline

Returns the size of the multibody system state vector x = [q; v] for this MultibodyPlant.

This will equal the number of generalized positions (see num_positions()) plus the number of generalized velocities (see num_velocities()). Notice however that the state of a MultibodyPlant, stored in its Context, can actually contain other variables such as integrated power and discrete states.

◆ num_positions() [1/2]

int num_positions ( ) const
inline

Returns the size of the generalized position vector q for this MultibodyPlant.

◆ num_positions() [2/2]

int num_positions ( ModelInstanceIndex  model_instance) const
inline

Returns the size of the generalized position vector q for a specific model instance.

◆ num_velocities() [1/2]

int num_velocities ( ) const
inline

Returns the size of the generalized velocity vector v for this MultibodyPlant.

◆ num_velocities() [2/2]

int num_velocities ( ModelInstanceIndex  model_instance) const
inline

Returns the size of the generalized velocity vector v for a specific model instance.

◆ num_visual_geometries()

int num_visual_geometries ( ) const
inline

Returns the number of geometries registered for visualization.

This method can be called at any time during the lifetime of this plant, either pre- or post-finalize, see Finalize(). Post-finalize calls will always return the same value.

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ RegisterAsSourceForSceneGraph()

geometry::SourceId RegisterAsSourceForSceneGraph ( geometry::SceneGraph< T > *  scene_graph)

Registers this plant to serve as a source for an instance of SceneGraph.

This registration allows MultibodyPlant to register geometry with scene_graph for visualization and/or collision queries. Successive registration calls with SceneGraph must be performed on the same instance to which the pointer argument scene_graph points to. Failure to do so will result in runtime exceptions.

Parameters
scene_graphA valid non nullptr to the SceneGraph instance for which this plant will sever as a source, see SceneGraph documentation for further details.
Returns
the SourceId of this plant in scene_graph. It can also later on be retrieved with get_source_id().
Exceptions
std::exceptionif called post-finalize.
std::exceptionif scene_graph is the nullptr.
std::exceptionif called more than once.

◆ RegisterCollisionGeometry()

geometry::GeometryId RegisterCollisionGeometry ( const Body< T > &  body,
const Isometry3< double > &  X_BG,
const geometry::Shape shape,
const std::string &  name,
const CoulombFriction< double > &  coulomb_friction,
geometry::SceneGraph< T > *  scene_graph 
)

Registers geometry in a SceneGraph with a given geometry::Shape to be used for the contact modeling of a given body.

More than one geometry can be registered with a body, in which case the body's contact geometry is the union of all geometries registered to that body.

Parameters
[in]bodyThe body for which geometry is being registered.
[in]X_BGThe fixed pose of the geometry frame G in the body frame B.
[in]shapeThe geometry::Shape used for visualization. E.g.: geometry::Sphere, geometry::Cylinder, etc.
[in]coulomb_frictionCoulomb's law of friction coefficients to model friction on the surface of shape for the given body.
[out]scene_graphA valid, non-null pointer to a SceneGraph on which geometry will get registered.
Exceptions
std::exceptionif scene_graph is the nullptr.
std::exceptionif called post-finalize.
std::exceptionif scene_graph does not correspond to the same instance with which RegisterAsSourceForSceneGraph() was called.

◆ RegisterVisualGeometry() [1/2]

geometry::GeometryId RegisterVisualGeometry ( const Body< T > &  body,
const Isometry3< double > &  X_BG,
const geometry::Shape shape,
const std::string &  name,
const geometry::VisualMaterial material,
geometry::SceneGraph< T > *  scene_graph 
)

Registers geometry in a SceneGraph with a given geometry::Shape to be used for visualization of a given body.

Parameters
[in]bodyThe body for which geometry is being registered.
[in]X_BGThe fixed pose of the geometry frame G in the body frame B.
[in]shapeThe geometry::Shape used for visualization. E.g.: geometry::Sphere, geometry::Cylinder, etc.
[in]nameThe name for the geometry. It must satsify the requirements defined in drake::geometry::GeometryInstance.
[in]materialThe visual material to assign to the geometry.
[out]scene_graphA valid non nullptr to a SceneGraph on which geometry will get registered.
Exceptions
std::exceptionif scene_graph is the nullptr.
std::exceptionif called post-finalize.
std::exceptionif scene_graph does not correspond to the same instance with which RegisterAsSourceForSceneGraph() was called.
Returns
the id for the registered geometry.

◆ RegisterVisualGeometry() [2/2]

geometry::GeometryId RegisterVisualGeometry ( const Body< T > &  body,
const Isometry3< double > &  X_BG,
const geometry::Shape shape,
const std::string &  name,
geometry::SceneGraph< T > *  scene_graph 
)

Overload for visual geometry registration; it implicitly assigns the default material.

◆ set_penetration_allowance()

void set_penetration_allowance ( double  penetration_allowance = 0.001)

Sets the penetration allowance used to estimate the coefficients in the penalty method used to impose non-penetration among bodies.

Refer to the section Contact by penalty method for further details.

◆ set_stiction_tolerance()

void set_stiction_tolerance ( double  v_stiction = 0.001)
inline

Sets the stiction tolerance v_stiction for the Stribeck model, where v_stiction must be specified in m/s (meters per second.) v_stiction defaults to a value of 1 millimeter per second.

Exceptions
std::exceptionif v_stiction is non-positive.

◆ SetDefaultState()

void SetDefaultState ( const systems::Context< T > &  context,
systems::State< T > *  state 
) const
inlineoverridevirtual

Sets the state in context so that generalized positions and velocities are zero.

Exceptions
std::exceptionif called pre-finalize. See Finalize().

Implements System< T >.

Reimplemented in AxiallySymmetricFreeBodyPlant< T >.

◆ time_step()

double time_step ( ) const
inline

The time step (or period) used to model this plant as a discrete system with periodic updates.

Returns 0 (zero) if the plant is modeled as a continuous system. This property of the plant is specified at construction and therefore this query can be performed either pre- or post- finalize, see Finalize().

See also
MultibodyPlant::MultibodyPlant(double)

◆ WeldFrames()

const WeldJoint< T > & WeldFrames ( const Frame< T > &  A,
const Frame< T > &  B,
const Isometry3< double > &  X_AB = Isometry3<double>::Identity() 
)

Welds frames A and B with relative pose X_AB.

That is, the pose of frame B in frame A is fixed, with value X_AB. The call to this method creates and adds a new WeldJoint to the model. The new WeldJoint is named as: A.name() + "_welds_to_" + B.name().

Returns
a constant reference to the WeldJoint welding frames A and B.

◆ world_body()

const RigidBody<T>& world_body ( ) const
inline

Returns a constant reference to the world body.

◆ world_frame()

const BodyFrame<T>& world_frame ( ) const
inline

Returns a constant reference to the world frame.

Friends And Related Function Documentation

◆ MultibodyPlant

friend class MultibodyPlant
friend

◆ MultibodyPlantTester

friend class MultibodyPlantTester
friend

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