Drake
Drake C++ Documentation
DeformableModel< T > Class Template Referencefinal

Detailed Description

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

DeformableModel implements the interface in PhysicalModel and provides the functionalities to specify deformable bodies.

Unlike rigid bodies, the shape of deformable bodies can change in a simulation. Each deformable body is modeled as a volumetric mesh with persisting topology, changing vertex positions, and an approximated signed distance field. A finite element model is built for each registered deformable body that is used to evaluate the dynamics of the body.

Warning
This feature is considered to be experimental and may change or be removed at any time, without any deprecation notice ahead of time.
Template Parameters
TThe scalar type, which must be double.

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

Public Member Functions

 DeformableModel (MultibodyPlant< T > *plant)
 (Internal only) Constructs a DeformableModel to be owned by the given MultibodyPlant. More...
 
 ~DeformableModel () final
 
int num_bodies () const
 Returns the number of deformable bodies registered with this DeformableModel. More...
 
DeformableBodyId RegisterDeformableBody (std::unique_ptr< geometry::GeometryInstance > geometry_instance, ModelInstanceIndex model_instance, const fem::DeformableBodyConfig< T > &config, double resolution_hint)
 Registers a deformable body in this DeformableModel with the given GeometryInstance. More...
 
DeformableBodyId RegisterDeformableBody (std::unique_ptr< geometry::GeometryInstance > geometry_instance, const fem::DeformableBodyConfig< T > &config, double resolution_hint)
 Registers a deformable body in this DeformableModel with the default model instance. More...
 
void SetWallBoundaryCondition (DeformableBodyId id, const Vector3< T > &p_WQ, const Vector3< T > &n_W)
 Sets wall boundary conditions for the body with the given id. More...
 
MultibodyConstraintId AddFixedConstraint (DeformableBodyId body_A_id, const RigidBody< T > &body_B, const math::RigidTransform< double > &X_BA, const geometry::Shape &shape_G, const math::RigidTransform< double > &X_BG)
 Defines a fixed constraint between a deformable body A and a rigid body B. More...
 
systems::DiscreteStateIndex GetDiscreteStateIndex (DeformableBodyId id) const
 Returns the discrete state index of the deformable body identified by the given id. More...
 
void SetPositions (systems::Context< T > *context, DeformableBodyId id, const Eigen::Ref< const Matrix3X< T >> &q) const
 Sets the vertex positions of the deformable body with the given id in the provided context. More...
 
void SetVelocities (systems::Context< T > *context, DeformableBodyId id, const Eigen::Ref< const Matrix3X< T >> &v) const
 Sets the vertex velocities of the deformable body with the given id in the provided context. More...
 
void SetPositionsAndVelocities (systems::Context< T > *context, DeformableBodyId id, const Eigen::Ref< const Matrix3X< T >> &q, const Eigen::Ref< const Matrix3X< T >> &v) const
 Sets the vertex positions and velocities of the deformable body with the given id in the provided context. More...
 
Matrix3X< T > GetPositions (const systems::Context< T > &context, DeformableBodyId id) const
 Returns the matrix of vertex positions for the deformable body with the given id in the provided context. More...
 
Matrix3X< T > GetVelocities (const systems::Context< T > &context, DeformableBodyId id) const
 Returns the matrix of vertex velocities for the deformable body with the given id in the provided context. More...
 
Matrix3X< T > GetPositionsAndVelocities (const systems::Context< T > &context, DeformableBodyId id) const
 Returns the matrix of vertex positions and velocities for the deformable body with the given id in the provided context. More...
 
void AddExternalForce (std::unique_ptr< ForceDensityFieldBase< T >> external_force)
 Registers an external force density field that applies external force to all deformable bodies. More...
 
const std::vector< const ForceDensityFieldBase< T > * > & GetExternalForces (DeformableBodyId id) const
 Returns the force density fields acting on the deformable body with the given id. More...
 
void Disable (DeformableBodyId id, systems::Context< T > *context) const
 Disables the deformable body with the given id in the given context. More...
 
void Enable (DeformableBodyId id, systems::Context< T > *context) const
 Enables the deformable body with the given id in the given context. More...
 
bool is_enabled (DeformableBodyId id, const systems::Context< T > &context) const
 
const fem::FemModel< T > & GetFemModel (DeformableBodyId id) const
 Returns the FemModel for the body with id. More...
 
const VectorX< double > & GetReferencePositions (DeformableBodyId id) const
 Returns the reference positions of the vertices of the deformable body identified by the given id. More...
 
DeformableBodyId GetBodyId (DeformableBodyIndex index) const
 Returns the DeformableBodyId of the body with the given body index. More...
 
DeformableBodyId GetBodyId (geometry::GeometryId geometry_id) const
 Returns the DeformableBodyId associated with the given geometry_id. More...
 
const DeformableBody< T > & GetBody (DeformableBodyId id) const
 Returns the deformable body with the given id. More...
 
const DeformableBody< T > & GetBody (DeformableBodyIndex index) const
 Returns the deformable body with the given index. More...
 
DeformableBody< T > & GetMutableBody (DeformableBodyId id)
 Returns a mutable reference to the deformable body with the given id. More...
 
bool HasBodyNamed (const std::string &name) const
 Returns true if and only if a deformable body with the given name has been registered with this model. More...
 
bool HasBodyNamed (const std::string &name, ModelInstanceIndex model_instance) const
 Returns true if and only if a deformable body with the given name has been registered with this model under the given model_instance. More...
 
const DeformableBody< T > & GetBodyByName (const std::string &name) const
 Returns the DeformableBody with the given name. More...
 
const DeformableBody< T > & GetBodyByName (const std::string &name, ModelInstanceIndex model_instance) const
 Returns the DeformableBody with the given name from the given model instance. More...
 
std::vector< DeformableBodyIdGetBodyIds (ModelInstanceIndex model_instance) const
 Returns the DeformableIds of the bodies that belong to the given model instance. More...
 
DeformableBodyIndex GetBodyIndex (DeformableBodyId id) const
 Returns the DeformableBodyIndex of the body with the given id. More...
 
geometry::GeometryId GetGeometryId (DeformableBodyId id) const
 Returns the GeometryId of the geometry associated with the body with the given id. More...
 
bool HasConstraint (DeformableBodyId id) const
 Returns the true if the deformable body with the given id has constraints associated with it. More...
 
const multibody::fem::internal::DiscreteTimeIntegrator< T > & integrator () const
 (Internal use only) Returns the time integrator used to for all FemModels in this model. More...
 
systems::OutputPortIndex configuration_output_port_index () const
 (Internal use only) Returns the output port index of the vertex positions port for all registered deformable bodies. More...
 
bool is_empty () const
 Returns true if there's no deformable body or external force registered to this DeformableModel. More...
 
bool is_cloneable_to_double () const final
 
bool is_cloneable_to_autodiff () const final
 Returns true if and only if this DeformableModel is empty. More...
 
bool is_cloneable_to_symbolic () const final
 Returns true if and only if this DeformableModel is empty. More...
 
void SetParallelism (Parallelism parallelism)
 (Internal use only) Configures the parallelism that this DeformableModel uses when opportunities for parallel computation arises. More...
 
Parallelism parallelism () const
 (Internal use only) Returns the parallelism that this DeformableModel uses when opportunities for parallel computation arises. More...
 
void SetDefaultState (const systems::Context< T > &context, systems::State< T > *state) const
 (Internal use only) Sets the default state for the deformable model. More...
 
Does not allow copy, move, or assignment
 DeformableModel (const DeformableModel &)=delete
 
DeformableModeloperator= (const DeformableModel &)=delete
 
 DeformableModel (DeformableModel &&)=delete
 
DeformableModeloperator= (DeformableModel &&)=delete
 
- Public Member Functions inherited from PhysicalModel< T >
 PhysicalModel (MultibodyPlant< T > *owning_plant)
 Constructs a PhysicalModel owned by the given owning_plant. More...
 
 ~PhysicalModel () override
 
template<typename ScalarType >
std::unique_ptr< PhysicalModel< ScalarType > > CloneToScalar (MultibodyPlant< ScalarType > *plant) const
 (Internal only) Creates a clone of this concrete PhysicalModel object with the scalar type ScalarType to be owned by the given plant. More...
 
bool is_cloneable_to_double () const override
 Defaults to false. More...
 
bool is_cloneable_to_autodiff () const override
 Defaults to false. More...
 
bool is_cloneable_to_symbolic () const override
 Defaults to false. More...
 
void DeclareSystemResources ()
 (Internal only) MultibodyPlant calls this from within Finalize() to declare additional system resources. More...
 
void DeclareSceneGraphPorts ()
 (Internal only) Declares zero or more output ports in the owning MultibodyPlant to communicate with a SceneGraph. More...
 
PhysicalModelPointerVariant< T > ToPhysicalModelPointerVariant () const
 Returns (a const pointer to) the specific model variant of this PhysicalModel. More...
 
const MultibodyPlant< T > & plant () const
 
 PhysicalModel (const PhysicalModel &)=delete
 
PhysicalModeloperator= (const PhysicalModel &)=delete
 
 PhysicalModel (PhysicalModel &&)=delete
 
PhysicalModeloperator= (PhysicalModel &&)=delete
 

Friends

template<typename U >
class DeformableModel
 

Additional Inherited Members

- Protected Member Functions inherited from PhysicalModel< T >
MultibodyPlant< T > * mutable_plant ()
 
const internal::MultibodyTree< T > & internal_tree () const
 
void ThrowIfSystemResourcesDeclared (const char *function_name) const
 
void ThrowIfSystemResourcesNotDeclared (const char *function_name) const
 
geometry::SceneGraph< T > & mutable_scene_graph ()
 
systems::DiscreteStateIndex DeclareDiscreteState (const VectorX< T > &model_value)
 
systems::AbstractParameterIndex DeclareAbstractParameter (const AbstractValue &model_value)
 
systems::LeafOutputPort< T > & DeclareAbstractOutputPort (std::string name, typename systems::LeafOutputPort< T >::AllocCallback alloc_function, typename systems::LeafOutputPort< T >::CalcCallback calc_function, std::set< systems::DependencyTicket > prerequisites_of_calc={ systems::System< T >::all_sources_ticket()})
 
systems::LeafOutputPort< T > & DeclareVectorOutputPort (std::string name, const systems::BasicVector< T > &model_vector, typename systems::LeafOutputPort< T >::CalcVectorCallback vector_calc_function, std::set< systems::DependencyTicket > prerequisites_of_calc={ systems::System< T >::all_sources_ticket()})
 

Constructor & Destructor Documentation

◆ DeformableModel() [1/3]

DeformableModel ( const DeformableModel< T > &  )
delete

◆ DeformableModel() [2/3]

DeformableModel ( DeformableModel< T > &&  )
delete

◆ DeformableModel() [3/3]

DeformableModel ( MultibodyPlant< T > *  plant)
explicit

(Internal only) Constructs a DeformableModel to be owned by the given MultibodyPlant.

This constructor is only intended to be called internally by MultibodyPlant.

Precondition
plant != nullptr.
Finalize() has not been called on plant.

◆ ~DeformableModel()

~DeformableModel ( )
final

Member Function Documentation

◆ AddExternalForce()

void AddExternalForce ( std::unique_ptr< ForceDensityFieldBase< T >>  external_force)

Registers an external force density field that applies external force to all deformable bodies.

Exceptions
std::exceptionif this DeformableModel is not of scalar type double.
std::exceptionif this DeformableModel belongs to a continuous MultibodyPlant.
std::exceptionif Finalize() has been called on the multibody plant owning this deformable model.

◆ AddFixedConstraint()

MultibodyConstraintId AddFixedConstraint ( DeformableBodyId  body_A_id,
const RigidBody< T > &  body_B,
const math::RigidTransform< double > &  X_BA,
const geometry::Shape shape_G,
const math::RigidTransform< double > &  X_BG 
)

Defines a fixed constraint between a deformable body A and a rigid body B.

Such a fixed constraint is modeled as distance holonomic constraints:

p_PᵢQᵢ(q) = 0 for each constrained vertex Pᵢ

where Pᵢ is the i-th vertex of the deformable body under constraint and Qᵢ is a point rigidly affixed to the rigid body B. To specify the constraint, we put the reference mesh M of the deformable body A in B's body frame with the given pose X_BA and prescribe a shape G with pose X_BG in B's body frame. All vertices Pᵢ in M that are inside (or on the surface of) G are subject to the fixed constraints with Qᵢ being coincident with Pᵢ when M is in pose X_BA. p_PᵢQᵢ(q) denotes the relative position of point Qᵢ with respect to point Pᵢ as a function of the configuration of the model q. Imposing this constraint forces Pᵢ and Qᵢ to be coincident for each vertex i of the deformable body specified to be under constraint.

Parameters
[in]body_A_idThe unique id of the deformable body under the fixed constraint.
[in]body_BThe rigid body under constraint.
[in]X_BAThe pose of deformable body A's reference mesh in B's body frame
[in]shape_GThe prescribed geometry shape, attached to rigid body B, used to determine which vertices of the deformable body A is under constraint.
[in]X_BGThe fixed pose of the geometry frame of the given shape in body B's frame.
Returns
the unique id of the newly added constraint.
Exceptions
std::exceptionif no deformable body with the given body_A_id has been registered.
std::exceptionunless body_B is registered with the same multibody plant owning this deformable model.
std::exceptionif Finalize() has been called on the multibody plant owning this deformable model.
std::exceptionif this DeformableModel is not of scalar type double.
std::exceptionif no constraint is added (i.e. no vertex of the deformable body is inside the given shape with the given poses).

◆ configuration_output_port_index()

systems::OutputPortIndex configuration_output_port_index ( ) const

(Internal use only) Returns the output port index of the vertex positions port for all registered deformable bodies.

Exceptions
std::exceptionif called before DeclareSceneGraphPorts() is called.

◆ Disable()

void Disable ( DeformableBodyId  id,
systems::Context< T > *  context 
) const

Disables the deformable body with the given id in the given context.

Disabling a deformable body sets its vertex velocities and accelerations to zero and freezes its vertex positions. A disabled deformable body is not subject to any constraint (e.g. frictional contact constraint or fixed constraint); it does not move under the influence of external forces (e.g. gravity); and it does not necessarily satisfy the prescribed boundary condition (if any). On the flip side, a disabled deformable body does not affect the dynamics of other bodies, even if the collision between the disabled body's geometry and other geometries is not filtered. Effectively, the physics of the deformable body stop being computed. The deformable body can be enabled by calling Enable(). Calling Disable() on a body which is already disabled has no effect.

See also
Enable().
Exceptions
std::exceptionif the passed in context isn't compatible with the MultibodyPlant associated with this DeformableModel.
std::exceptionif a deformable body with the given id is not registered.
std::exceptionif context is null.

◆ Enable()

void Enable ( DeformableBodyId  id,
systems::Context< T > *  context 
) const

Enables the deformable body with the given id in the given context.

Calling Enable() on a body which is already enabled has no effect.

See also
Disable().
Exceptions
std::exceptionif the passed in context isn't compatible with the MultibodyPlant associated with this DeformableModel.
std::exceptionif a deformable body with the given id is not registered.
std::exceptionif context is null.

◆ GetBody() [1/2]

const DeformableBody<T>& GetBody ( DeformableBodyId  id) const

Returns the deformable body with the given id.

Exceptions
std::exceptionif no deformable body with the given id has been registered in this model.

◆ GetBody() [2/2]

const DeformableBody<T>& GetBody ( DeformableBodyIndex  index) const

Returns the deformable body with the given index.

Exceptions
std::exceptionif no deformable body with the given index is registered in this model.

◆ GetBodyByName() [1/2]

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

Returns the DeformableBody with the given name.

Exceptions
std::exceptionif there's no body with the given name or if more than one model instance contains a deformable body with the given name.

◆ GetBodyByName() [2/2]

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

Returns the DeformableBody with the given name from the given model instance.

Exceptions
std::exceptionif there's no body with the given name that is registered with the given model instance.

◆ GetBodyId() [1/2]

DeformableBodyId GetBodyId ( DeformableBodyIndex  index) const

Returns the DeformableBodyId of the body with the given body index.

Exceptions
std::exceptionif no deformable body with the given index has been registered in this model.

◆ GetBodyId() [2/2]

DeformableBodyId GetBodyId ( geometry::GeometryId  geometry_id) const

Returns the DeformableBodyId associated with the given geometry_id.

Exceptions
std::exceptionif the given geometry_id does not correspond to a deformable body registered with this model.

◆ GetBodyIds()

std::vector<DeformableBodyId> GetBodyIds ( ModelInstanceIndex  model_instance) const

Returns the DeformableIds of the bodies that belong to the given model instance.

Returns the empty vector if no deformable bodies are registered with the given model instance.

◆ GetBodyIndex()

DeformableBodyIndex GetBodyIndex ( DeformableBodyId  id) const

Returns the DeformableBodyIndex of the body with the given id.

Exceptions
std::exceptionif no body with the given id has been registered.

◆ GetDiscreteStateIndex()

systems::DiscreteStateIndex GetDiscreteStateIndex ( DeformableBodyId  id) const

Returns the discrete state index of the deformable body identified by the given id.

Exceptions
std::exceptionif MultibodyPlant::Finalize() has not been called yet. or if no deformable body with the given id has been registered in this model.

◆ GetExternalForces()

const std::vector<const ForceDensityFieldBase<T>*>& GetExternalForces ( DeformableBodyId  id) const

Returns the force density fields acting on the deformable body with the given id.

Exceptions
std::exceptionif MultibodyPlant::Finalize() has not been called yet. or if no deformable body with the given id has been registered in this model.

◆ GetFemModel()

const fem::FemModel<T>& GetFemModel ( DeformableBodyId  id) const

Returns the FemModel for the body with id.

Exceptions
exceptionif no deformable body with id is registered with this DeformableModel.

◆ GetGeometryId()

geometry::GeometryId GetGeometryId ( DeformableBodyId  id) const

Returns the GeometryId of the geometry associated with the body with the given id.

Exceptions
std::exceptionif no body with the given id has been registered.

◆ GetMutableBody()

DeformableBody<T>& GetMutableBody ( DeformableBodyId  id)

Returns a mutable reference to the deformable body with the given id.

Exceptions
std::exceptionif no deformable body with the given id has been registered in this model.

◆ GetPositions()

Matrix3X<T> GetPositions ( const systems::Context< T > &  context,
DeformableBodyId  id 
) const

Returns the matrix of vertex positions for the deformable body with the given id in the provided context.

Parameters
[in]contextThe context associated with the MultibodyPlant that owns this DeformableModel.
[in]idThe identifier of the deformable body whose positions are being queried.
Return values
qA 3×N matrix containing the positions of all vertices of the body.
Exceptions
std::exceptionif any of the following conditions are met:
  1. context does not belong to the MultibodyPlant associated with this DeformableModel.
  2. No body with the given id is registered.
  3. Finalize() has not been called on the MultibodyPlant that owns this deformable model.

◆ GetPositionsAndVelocities()

Matrix3X<T> GetPositionsAndVelocities ( const systems::Context< T > &  context,
DeformableBodyId  id 
) const

Returns the matrix of vertex positions and velocities for the deformable body with the given id in the provided context.

The first N columns are the positions and the next N columns are the velocities.

Parameters
[in]contextThe context associated with the MultibodyPlant that owns this DeformableModel.
[in]idThe identifier of the deformable body whose state is being queried.
Returns
A 3x2N matrix containing the positions and velocities of all vertices of the body.
Exceptions
std::exceptionif any of the following conditions are met:
  1. context does not belong to the MultibodyPlant associated with this DeformableModel.
  2. No body with the given id is registered.
  3. Finalize() has not been called on the MultibodyPlant that owns this deformable model.

◆ GetReferencePositions()

const VectorX<double>& GetReferencePositions ( DeformableBodyId  id) const

Returns the reference positions of the vertices of the deformable body identified by the given id.

The reference positions are represented as a VectorX with 3N values where N is the number of vertices. The x-, y-, and z-positions (measured and expressed in the world frame) of the j-th vertex are 3j, 3j + 1, and 3j + 2 in the VectorX.

Exceptions
std::exceptionif no deformable body with the given id has been registered in this model.

◆ GetVelocities()

Matrix3X<T> GetVelocities ( const systems::Context< T > &  context,
DeformableBodyId  id 
) const

Returns the matrix of vertex velocities for the deformable body with the given id in the provided context.

Parameters
[in]contextThe context associated with the MultibodyPlant that owns this DeformableModel.
[in]idThe identifier of the deformable body whose velocities are being queried.
Return values
vA 3×N matrix containing the velocities of all vertices of the body.
Exceptions
std::exceptionif any of the following conditions are met:
  1. context does not belong to the MultibodyPlant associated with this DeformableModel.
  2. No body with the given id is registered.
  3. Finalize() has not been called on the MultibodyPlant that owns this deformable model.

◆ HasBodyNamed() [1/2]

bool HasBodyNamed ( const std::string &  name) const

Returns true if and only if a deformable body with the given name has been registered with this model.

◆ HasBodyNamed() [2/2]

bool HasBodyNamed ( const std::string &  name,
ModelInstanceIndex  model_instance 
) const

Returns true if and only if a deformable body with the given name has been registered with this model under the given model_instance.

◆ HasConstraint()

bool HasConstraint ( DeformableBodyId  id) const

Returns the true if the deformable body with the given id has constraints associated with it.

◆ integrator()

const multibody::fem::internal::DiscreteTimeIntegrator<T>& integrator ( ) const

(Internal use only) Returns the time integrator used to for all FemModels in this model.

Exceptions
std::exceptionif the integrator hasn't been set.

◆ is_cloneable_to_autodiff()

bool is_cloneable_to_autodiff ( ) const
final

Returns true if and only if this DeformableModel is empty.

◆ is_cloneable_to_double()

bool is_cloneable_to_double ( ) const
final

◆ is_cloneable_to_symbolic()

bool is_cloneable_to_symbolic ( ) const
final

Returns true if and only if this DeformableModel is empty.

◆ is_empty()

bool is_empty ( ) const

Returns true if there's no deformable body or external force registered to this DeformableModel.

◆ is_enabled()

bool is_enabled ( DeformableBodyId  id,
const systems::Context< T > &  context 
) const
Returns
true if and only if the deformable body with the given id is enabled.
Exceptions
std::exceptionif the passed in context isn't compatible with the MultibodyPlant associated with this DeformableModel.
std::exceptionif a deformable body with the given id is not registered.

◆ num_bodies()

int num_bodies ( ) const

Returns the number of deformable bodies registered with this DeformableModel.

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ parallelism()

Parallelism parallelism ( ) const

(Internal use only) Returns the parallelism that this DeformableModel uses when opportunities for parallel computation arises.

◆ RegisterDeformableBody() [1/2]

DeformableBodyId RegisterDeformableBody ( std::unique_ptr< geometry::GeometryInstance geometry_instance,
ModelInstanceIndex  model_instance,
const fem::DeformableBodyConfig< T > &  config,
double  resolution_hint 
)

Registers a deformable body in this DeformableModel with the given GeometryInstance.

The body is represented in the world frame and simulated with FEM with linear elements and a first order quadrature rule that integrates linear functions exactly. See FemModel for details. Returns a unique identifier for the added geometry.

Parameters
[in]geometry_instanceThe geometry to be registered with the model.
[in]configThe physical properties of deformable body.
[in]model_instanceThe model instance index which this body is part of.
[in]resolution_hintThe parameter that guides the level of mesh refinement of the deformable geometry. It has length units (in meters) and roughly corresponds to a typical edge length in the resulting mesh for a primitive shape.
Precondition
resolution_hint > 0.
Exceptions
std::exceptionif this DeformableModel is not of scalar type double.
std::exceptionif this DeformableModel belongs to a continuous MultibodyPlant.
std::exceptionif the model instance does not exist.
std::exceptionif a deformable body with the same name has already been registered to the model instance.
std::exceptionif Finalize() has been called on the multibody plant owning this deformable model.

◆ RegisterDeformableBody() [2/2]

DeformableBodyId RegisterDeformableBody ( std::unique_ptr< geometry::GeometryInstance geometry_instance,
const fem::DeformableBodyConfig< T > &  config,
double  resolution_hint 
)

Registers a deformable body in this DeformableModel with the default model instance.

◆ SetDefaultState()

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

(Internal use only) Sets the default state for the deformable model.

This should only be called by MultibodyPlant as a part of MultibodyPlant::SetDefaultState().

◆ SetParallelism()

void SetParallelism ( Parallelism  parallelism)

(Internal use only) Configures the parallelism that this DeformableModel uses when opportunities for parallel computation arises.

◆ SetPositions()

void SetPositions ( systems::Context< T > *  context,
DeformableBodyId  id,
const Eigen::Ref< const Matrix3X< T >> &  q 
) const

Sets the vertex positions of the deformable body with the given id in the provided context.

Parameters
[in,out]contextThe context associated with the MultibodyPlant that owns this DeformableModel.
[in]idThe identifier of the deformable body whose positions are being set.
[in]qA 3×N matrix of vertex positions.
Exceptions
std::exceptionif any of the following conditions are met:
  1. context is nullptr.
  2. context does not belong to the MultibodyPlant associated with this DeformableModel.
  3. No body with the given id is registered.
  4. The number of columns of q does not match the number of vertices of the body.
  5. q contains non-finite values.
  6. Finalize() has not been called on the MultibodyPlant that owns this deformable model.

◆ SetPositionsAndVelocities()

void SetPositionsAndVelocities ( systems::Context< T > *  context,
DeformableBodyId  id,
const Eigen::Ref< const Matrix3X< T >> &  q,
const Eigen::Ref< const Matrix3X< T >> &  v 
) const

Sets the vertex positions and velocities of the deformable body with the given id in the provided context.

Parameters
[in,out]contextThe context associated with the MultibodyPlant that owns this DeformableModel.
[in]idThe identifier of the deformable body whose positions and velocities are being set.
[in]qA 3×N matrix of vertex positions.
[in]vA 3×N matrix of vertex velocities.
Exceptions
std::exceptionif any of the following conditions are met:
  1. context is nullptr.
  2. context does not belong to the MultibodyPlant associated with this DeformableModel.
  3. No body with the given id is registered.
  4. The number of columns of q or v does not match the number of vertices of the body.
  5. q or v contains non-finite values.
  6. Finalize() has not been called on the MultibodyPlant that owns this deformable model.

◆ SetVelocities()

void SetVelocities ( systems::Context< T > *  context,
DeformableBodyId  id,
const Eigen::Ref< const Matrix3X< T >> &  v 
) const

Sets the vertex velocities of the deformable body with the given id in the provided context.

Parameters
[in,out]contextThe context associated with the MultibodyPlant that owns this DeformableModel.
[in]idThe identifier of the deformable body whose velocities are being set.
[in]vA 3×N matrix of vertex velocities.
Exceptions
std::exceptionif any of the following conditions are met:
  1. context is nullptr.
  2. context does not belong to the MultibodyPlant associated with this DeformableModel.
  3. No body with the given id is registered.
  4. The number of columns of v does not match the number of vertices of the body.
  5. v contains non-finite values.
  6. Finalize() has not been called on the MultibodyPlant that owns this deformable model.

◆ SetWallBoundaryCondition()

void SetWallBoundaryCondition ( DeformableBodyId  id,
const Vector3< T > &  p_WQ,
const Vector3< T > &  n_W 
)

Sets wall boundary conditions for the body with the given id.

All vertices of the mesh of the deformable body whose reference positions are inside the prescribed open half space are put under zero displacement boundary conditions. The open half space is defined by a plane with outward normal n_W. A vertex V is considered to be subject to the boundary condition if n̂ ⋅ p_QV < 0 where Q is a point on the plane and n̂ is normalized n_W.

Parameters
[in]idThe body to be put under boundary condition.
[in]p_WQThe position of a point Q on the plane in the world frame.
[in]n_WOutward normal to the half space expressed in the world frame.
Precondition
n_W.norm() > 1e-10.
Warning
Be aware of round-off errors in floating computations when placing a vertex very close to the plane defining the half space.
Exceptions
std::exceptionif no deformable body with the given id has been registered in this model.

Friends And Related Function Documentation

◆ DeformableModel

friend class DeformableModel
friend

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