Go to the documentation of this file. 1 #ifndef SimTK_SIMBODY_MOBILIZED_BODY_CUSTOM_H_
2 #define SimTK_SIMBODY_MOBILIZED_BODY_CUSTOM_H_
75 class ImplementationImpl;
119 #ifndef SimTK_SIMBODY_DEFINING_MOBILIZED_BODY
121 MobilizedBody::Custom::ImplementationImpl>;
131 :
public PIMPLHandle<Implementation,ImplementationImpl>
300 int nIn,
const Real* in,
int nOut,
Real* out)
const;
313 int nIn,
const Real* in,
int nOut,
Real* out)
const;
326 int nIn,
const Real* in,
int nOut,
Real* out)
const;
384 virtual void realizeTopology(
State&)
const { }
455 friend class MobilizedBody::CustomImpl;
460 #endif // SimTK_SIMBODY_MOBILIZED_BODY_CUSTOM_H_
This object is intended to contain all state information for a SimTK::System, except topological info...
Definition: State.h:280
virtual void realizeDynamics(const State &) const
The Matter Subsystem's realizeDynamics() method will call this method along with the built-in Mobiliz...
Definition: MobilizedBody_Custom.h:437
Vector getQDot(const State &s) const
Return a Vector containing all the generalized coordinate derivatives qdot currently in use by this m...
Custom(MobilizedBody &parent, Implementation *implementation, const Transform &X_PF, const Body &bodyInfo, const Transform &X_BM, Direction direction=Forward)
Create a Custom mobilizer between an existing parent (inboard) body P and a new child (outboard) body...
#define SimTK_INSERT_DERIVED_HANDLE_DECLARATIONS(DERIVED, DERIVED_IMPL, PARENT)
Definition: PrivateImplementation.h:343
Vector getQDotDot(const State &s) const
Return a Vector containing all the generalized coordinate second derivatives qdotdot currently in use...
virtual void multiplyByNDot(const State &s, bool transposeMatrix, int nIn, const Real *in, int nOut, Real *out) const
Calculate out_q = NDot(q)*in_u or out_u = ~NDot*in_q.
virtual void multiplyByHDotTranspose(const State &s, const SpatialVec &F, int nu, Real *f) const =0
Calculate f = ~HDot*F where F is a spatial vector and f is its mapping onto the mobilities.
virtual void realizeVelocity(const State &) const
The Matter Subsystem's realizeVelocity() method will call this method along with the built-in Mobiliz...
Definition: MobilizedBody_Custom.h:429
This is a System that represents the dynamics of a particle moving along a smooth surface.
Definition: Assembler.h:37
Transform getMobilizerTransform(const State &s) const
Get the cross-mobilizer transform X_FM, the body's "moving" mobilizer frame M measured and expressed ...
The handle class MobilizedBody::Custom (dataless) and its companion class MobilizedBody::Custom::Impl...
Definition: MobilizedBody_Custom.h:72
virtual Transform calcMobilizerTransformFromQ(const State &s, int nq, const Real *q) const =0
Given values for this mobilizer's nq generalized coordinates q, compute X_FM(q), that is,...
virtual void realizeAcceleration(const State &) const
The Matter Subsystem's realizeAcceleration() method will call this method along with the built-in Mob...
Definition: MobilizedBody_Custom.h:445
This subsystem contains the bodies ("matter") in the multibody system, the mobilizers (joints) that d...
Definition: SimbodyMatterSubsystem.h:133
virtual void multiplyByHTranspose(const State &s, const SpatialVec &F, int nu, Real *f) const =0
Calculate f = ~H*F where F is a spatial force (torque+force) and f is its mapping onto the mobilities...
virtual SpatialVec multiplyByHDotMatrix(const State &s, int nu, const Real *u) const =0
Calculate A0_FM = HDot*u where HDot=HDot(q,u) is the time derivative of H.
virtual void realizeTime(const State &) const
The Matter Subsystem's realizeTime() method will call this method along with the built-in MobilizedBo...
Definition: MobilizedBody_Custom.h:413
virtual ~Implementation()
Destructor is virtual so derived classes get a chance to clean up if necessary.
Custom(MobilizedBody &parent, Implementation *implementation, const Body &bodyInfo, Direction direction=Forward)
Abbreviated constructor you can use if the mobilizer frames are coincident with the parent and child ...
virtual void calcDecorativeGeometryAndAppend(const State &s, Stage stage, Array_< DecorativeGeometry > &geom) const
Implement this optional method if you would like your MobilizedBody to generate any suggestions for g...
Definition: MobilizedBody_Custom.h:368
virtual void setQToFitTransform(const State &, const Transform &X_FM, int nq, Real *q) const
Find a set of q's for this mobilizer that best approximate the supplied Transform which requests a pa...
This class is basically a glorified enumerated type, type-safe and range checked but permitting conve...
Definition: Stage.h:66
virtual void realizeReport(const State &) const
The Matter Subsystem's realizeReport() method will call this method along with the built-in Mobilized...
Definition: MobilizedBody_Custom.h:452
virtual SpatialVec multiplyByHMatrix(const State &s, int nu, const Real *u) const =0
Calculate V_FM(u) = H*u where H=H(q) is the joint transition matrix mapping the mobilities to the rel...
void invalidateTopologyCache() const
Call this if you want to make sure that the next realizeTopology() call does something.
virtual void realizeModel(State &) const
The Matter Subsystem's realizeModel() method will call this method along with the built-in MobilizedB...
Definition: MobilizedBody_Custom.h:402
SpatialVec getMobilizerVelocity(const State &s) const
Get the cross-mobilizer velocity V_FM, the relative velocity of this body's "moving" mobilizer frame ...
This class provides some infrastructure useful in making SimTK Private Implementation (PIMPL) classes...
Definition: PrivateImplementation.h:106
A MobilizedBody is Simbody's fundamental body-and-joint object used to parameterize a system's motion...
Definition: MobilizedBody.h:169
Vector getU(const State &s) const
Return a Vector containing all the generalized speeds u currently in use by this mobilizer.
Vector getQ(const State &s) const
Return a Vector containing all the generalized coordinates q currently in use by this mobilizer.
Direction
Constructors can take an argument of this type to indicate that the mobilizer is being defined in the...
Definition: MobilizedBody.h:181
#define SimTK_SIMBODY_EXPORT
Definition: Simbody/include/simbody/internal/common.h:68
The Array_<T> container class is a plug-compatible replacement for the C++ standard template library ...
Definition: Array.h:53
SimTK_Real Real
This is the default compiled-in floating point type for SimTK, either float or double.
Definition: SimTKcommon/include/SimTKcommon/internal/common.h:606
bool getUseEulerAngles(const State &s) const
Get whether rotations are being represented as quaternions or Euler angles.
This is the implementation class for Custom mobilizers.
Definition: MobilizedBody_Custom.h:132
Implementation(SimbodyMatterSubsystem &, int nu, int nq, int nAngles=0)
This Implementation base class constructor sets the topological defaults for the number of mobilities...
virtual void realizePosition(const State &) const
The Matter Subsystem's realizePosition() method will call this method along with the built-in Mobiliz...
Definition: MobilizedBody_Custom.h:421
The Body class represents a reference frame that can be used to describe mass properties and geometry...
Definition: Body.h:55
virtual void multiplyByN(const State &s, bool transposeMatrix, int nIn, const Real *in, int nOut, Real *out) const
Calculate out_q = N(q)*in_u (e.g., qdot=N*u) or out_u = ~N*in_q.
Vector getUDot(const State &s) const
Return a Vector containing all the generalized accelerations udot currently in use by this mobilizer.
virtual void setUToFitVelocity(const State &, const SpatialVec &V_FM, int nu, Real *u) const
Find a set of u's (generalized speeds) for this mobilizer that best approximate the supplied spatial ...
virtual void realizeInstance(const State &) const
The Matter Subsystem's realizeInstance() method will call this method along with the built-in Mobiliz...
Definition: MobilizedBody_Custom.h:407
Custom()
Default constructor provides an empty handle that can be assigned to reference any MobilizedBody::Cus...
Definition: MobilizedBody_Custom.h:79
virtual void multiplyByNInv(const State &s, bool transposeMatrix, int nIn, const Real *in, int nOut, Real *out) const
Calculate out_u = NInv(q)*in_q (e.g., u=NInv*qdot) or out_q = ~NInv*in_u.
const Implementation & getImplementation() const
virtual Implementation * clone() const =0
This method should produce a deep copy identical to the concrete derived Implementation object underl...
Implementation & updImplementation()