Class CompoundBody

Contents

Class CompoundBody#

Nested Relationships#

Nested Types#

Inheritance Relationships#

Base Types#

Derived Type#

Class Documentation#

class CompoundBody : public Karana::Core::LockingBase, public Karana::Dynamics::BodyBase#

Represents a compound body.

This class is for compound bodies representing a subtree of bodies.

Subclassed by Karana::Dynamics::CECompoundBody

Public Functions

CompoundBody(std::string_view name, kc::ks_ptr<SubTree> parent_subtree, kc::ks_ptr<SubTree> component_bodies)#

CompoundBody constructor. The constructor is not meant to be called directly. Please use the create(…) method instead to create an instance.

Parameters:
  • name – body name

  • parent_subtree – the parent SubTree for the compound body

  • component_bodiesSubTree with the bodies behind aggregated

virtual ~CompoundBody()#

CompoundBody destructor.

inline virtual std::string_view typeString() const noexcept override#

Returns the type string of Base.

Returns:

The type string.

inline virtual const kc::id_t &id() const override#

Helper method for reconciling multiple inheritance.

Returns:

the object id

inline virtual std::string_view name() const override#

Helper method for reconciling multiple inheritance.

Returns:

the object name

virtual kc::ks_ptr<kf::FrameToFrame> f2f() const override#

Return the frame to frame for the CoordBase.

The returned f2f helps locate the CoordBase in the frames tree. Its transform etc data does not however represent the relative transformation changes from the changes to this object’s coordinates.

Returns:

Return the oframe Frame instance

inline const kc::ks_ptr<SubTree> &bodiesTree() const#

Return the SubTree for the aggregated bodies.

Returns:

the aggregated bodies SubTree

inline kc::ks_ptr<SubTree> physicalBodiesTree() const#

Return the SubTree with aggregated physical bodies.

The returned tree only has physical bodies, and is obtained by recursively expanding out any compound bodies within the aggregated bodies in this compound body

Returns:

the SubTree with all aggregated physical bodies

virtual void setQ(const Eigen::Ref<const km::Vec> &val) override#

Set the global (chart) Q coordinates.

Get/set methods for buffers for the generalized coordinates for this coordinate provider.

See Generalized Q, U etc coordinates section for discussion on Q, U etc coordinates.

Parameters:

val – Array of values.

virtual const km::Vec &getQ() const override#

Return the global (chart) Q coordinates.

See Generalized Q, U etc coordinates section for discussion on Q, U etc coordinates.

Returns:

Array of values.

virtual void setU(const Eigen::Ref<const km::Vec> &val) override#

Set the U velocity coordinates.

See Generalized Q, U etc coordinates section for discussion on Q, U etc coordinates.

Parameters:

val – Array of values.

virtual const km::Vec &getU() const override#

Return the U velocity coordinates.

See Generalized Q, U etc coordinates section for discussion on Q, U etc coordinates.

Returns:

Array of values.

virtual void setUdot(const Eigen::Ref<const km::Vec> &val) override#

Set the Udot acceleration coordinates.

See Generalized Q, U etc coordinates section for discussion on Q, U etc coordinates.

Parameters:

val – Array of values.

virtual const km::Vec &getUdot() const override#

Return the Udot acceleration coordinates.

See Generalized Q, U etc coordinates section for discussion on Q, U etc coordinates.

Returns:

Array of values.

virtual kc::ks_ptr<HingeBase> parentHinge() const override#

Return the parent hinge for the body.

For a PhysicalBody, parent hinge is a PhysicalHinge, while for a CompoundBody the parent hinge is a CompoundHinge

Returns:

The parent HingeBase instance

inline virtual const kc::ks_ptr<PhysicalBody> &physicalParentBody() const override#

Return the PhysicalBody inboard parent body for the body.

For a PhysicalBody this is simply the inboard body connected to it via a PhysicalHinge. For a CompoundBody, this is the single PhysicalBody to which the base pnodes in the CompoundBody are connected to.

Returns:

the parent PhysicalBody body instance

bool containsBody(const BodyBase &bd) const#

Return true if the specified BodyBase body is embedded at some level within this compound body.

Parameters:

bd – the input body

Returns:

true if the body belongs to the subtree

virtual std::string dumpString(std::string_view prefix, const Karana::Core::Base::DumpOptions *options) const override#

Return information about the object.

Parameters:
  • prefix – A string to use as prefix for each output line

  • options – Struct with options to tailor the output

Returns:

String with the information about the object.

inline const kc::ks_ptr<kc::DataCache<bool>> atbiMatrixCache() const#

Return the ATBI matrix cache.

Returns:

the ATBI matrix cache for the body

inline const kc::ks_ptr<kc::DataCache<bool>> atbiFilterCache() const#

Return the ATBI filter cache.

Returns:

the ATBI filter cache for the body

inline const kc::ks_ptr<kc::DataCache<bool>> atbiSmootherCache() const#

Return the ATBI smoother cache.

Returns:

the ATBI smoother cache for the body

inline const kc::ks_ptr<kc::DataCache<bool>> upsilonMatrixCache() const#

Return the OSCM related Upsilon cache.

Returns:

the OSCM related Upsilon cache for the body

inline const InvDynVectors &inverseDynamicsVectors() const#

Return the updated inverse dynamics vectors.

Returns:

the inverse dynamics vectors

km::Vec getInterBodyForceTreeFwdDyn() const#

Return the stacked vector of interbody forces for the embedded physical bodies (at their onodes)

Returns:

the onode referenced interbody forces vector

Public Static Functions

static kc::ks_ptr<CompoundBody> create(std::string_view name, kc::ks_ptr<SubTree> parent_subtree, kc::ks_ptr<SubTree> component_bodies)#

Factory create method for a CompoundBody.

Parameters:
  • parent_subtree – the parent SubTree for the compound body

  • name – body name

  • component_bodiesSubTree with the bodies behind aggregated

Returns:

the CompoundBody instance

Protected Functions

inline const kc::ks_ptr<kc::DataCache<InvDynVectors>> _inverseDynamicsCache() const#

Return the inverse dynamics cache.

Returns:

the inverse dynamics cache for the body

virtual kc::ks_ptr<CompoundSubhinge> _createSubhinge(kc::ks_ptr<CompoundHinge> hinge)#

Create the compound subhinge.

Parameters:

hinge – the parent compound hinge

Returns:

the new compound subhinge

void _discard(kc::ks_ptr<kc::Base> &base) override#

Parameters:

base – - Base pointer to the CompoundBody to discard.

std::vector<kc::ks_ptr<PhysicalBody>> _physicalBodiesList() const#

Return the list of physical bodies in the compound body.

Compute and return the list of physical bodies in the compound body. This does not include the physical root body

Returns:

list of physical bodies

km::Vec _getInterBodyForcePnodesTreeFwdDyn() const#

Return the stacked vector of interbody forces for the embedded physical bodies (at their pnodes)

Returns:

the pnode referenced interbody forces vector

virtual void _setupChildBodyCacheDependencies(kc::ks_ptr<BodyBase> child_body) override#

Set up cache dependencies with a child body.

Parameters:

child_body – the child body

virtual void _setupChildPhysicalBodyCacheDependencies(const PhysicalBody &child_body, bool skip_scatter) override#

Set up cache dependencies with a PhysicalBody child body.

Parameters:
  • child_body – the PhysicalBody child body

  • skip_scatter – if true, skip dependencies for the scatter phase

virtual void _setupChildCompoundBodyCacheDependencies(kc::ks_ptr<CompoundBody> child_body) override#

Set up cache dependencies with a CompoundBody child body.

Parameters:

child_body – the CompoundBody child body

virtual void _setupBaseBodyCacheDependencies() override#

Set up cache dependencies for a base body.

virtual void _teardownChildBodyCacheDependencies(kc::ks_ptr<BodyBase> child_body) override#

Tear down cache dependencies with a child body.

Parameters:

child_body – the child body

virtual void _teardownChildPhysicalBodyCacheDependencies(const PhysicalBody &child_body, bool skip_scatter) override#

Tear down cache dependencies with a PhysicalBody child body.

Parameters:
  • child_body – the PhysicalBody child body

  • skip_scatter – if true, skip dependencies for the scatter phase

virtual void _teardownChildCompoundBodyCacheDependencies(kc::ks_ptr<CompoundBody> child_body) override#

Tear down cache dependencies with a CompoundBody child body.

Parameters:

child_body – the CompoundBody child body

virtual void _teardownBaseBodyCacheDependencies() override#

Tear down cache dependencies for a base body.

void _initialSetup()#

Do initial setup of the compound body.

void _postSetup()#

Do post setup of the compound body.

void _oneTimeSetupDataCaches()#

One time setup of the body.

void _oneTimeTeardownDataCaches()#

One time tear down of the body.

km::Vec _compute_Phi_G_fprime(const km::Vec &f_prime) const#

compute Phi_G * f_prime (for the actual inter-body force)

Parameters:

f_prime – the input f’ forces

Returns:

the vector product

virtual km::Mat _oframe2pframePsi() const override#

Return the oframe to pframe Psi() matrix for this CoordBase.

For

  • a flex body this is a 6x(nU+6) psi matrix, the lhs is at the pnode, and the rhs at the body frame.

  • For a subhinge this is a 6x6 psi matrix, the lhs/rhs are the subhinge’s oframe and pframe pair.

  • a compound subhinge this is the 6x(nbodies*(6+nU)) psi matrix, the lhs is at the physical parent, and the rhs at the embedded bodies.

Returns:

the Psi() matrix

virtual km::Mat _oframe2pframePhi() const override#

Return the oframe to pframe Phi() matrix for this CoordBase.

For

  • a flex body this is a 6x(nU+6) phi matrix, the lhs is at the pnode, and the rhs at the body frame.

  • a subhinge this is a 6x6 phi matrix, the lhs/rhs are the subhinge’s oframe and pframe pair.

  • a compound subhinge this is the 6x(nbodies*(6+nU)) E_Phi_G matrix, the lhs is at the physical parent, and the rhs at the embedded bodies.

Returns:

the Phi() matrix

virtual km::Mat _pframe2otherPhi(const kf::Frame &other) const override#

Return the pframe to other frame Phi() matrix for this CoordBase.

For

  • a flex body this is a (nU+6)x6 phi matrix - the lhs is at the body frame, and the rhs at the other frame.

  • a subhinge this is a 6x6 matrix - the lhs is the pframe, and the rhs at the other frame.

  • a compound subhinge this is the (nbodies*(6+nU))x6 matrix, the lhs is at the embedded body frames, and the rhs at the other frame.

Parameters:

other – the other to frame

Returns:

the Phi() matrix

void _computeInvDynForces(InvDynVectors &val)#

Data cache callback for inter body force values for inverse dynamics.

Parameters:

valInvDynVectors with returned values.

virtual km::Vec _getFprimeTreeFwdDyn() const#

Compute the inter body force value for a base onode for inverse dynamics.

Compute fprime body forces using tree ATBI fwd dynamics data

Parameters:

base_onode – The hinge onode for the interbody forces

Returns:

the computed f’ stacked vector value

Protected Attributes

kc::ks_ptr<PhysicalBody> _physical_parent_body = nullptr#

the single physical virtual root parent body for the G subtree to use. Note that all the base bodies in the G subtree are attached (at possible multiple points) to the virtual root physical body even if the topological parent for this compound body may be a compound body itself. The aggregation principle, requires that all the base physical bodies in G physical have to be attached to a single physical body (at different onodes likely) and this parent physical body will belong to the parent compound body.

kc::RegistryList<HingeOnode> _base_onodes_list#

List of onodes for the bodies directly attached to the physical virtual root of the aggregated bodies. This list is used in the gather recursions step from the compound body into the parent body via the E_G operator.

kc::ks_ptr<SubTree> _bodies_tree = nullptr#

SubTree with the bodies in the compound body. This list may contain compound bodies as well.

kc::ks_ptr<SubTree> _physical_bodies_tree = nullptr#

SubTree with the physical bodies in the compound body. If _bodies_tree contains only physical bodies, then this variable simply points to it.

kc::ks_ptr<CompoundHinge> _hinge = nullptr#

the parent compound hinge for this body

std::map<kc::ks_ptr<PhysicalBody>, kc::RegistryList<HingeOnode>> _child_onodes_map#

Map of of child body attachment onodes for each aggregated body in the compound body. The list is the onodes for each of the aggregated bodies consists of just the subset of its child onodes that are connected to physical bodies external to the compound body.

std::map<kc::ks_ptr<PhysicalBody>, kc::ks_ptr<HingeOnode>> _base_onodes_map#

Map of of base onode for the branch each aggregated body in the compound body. All of the base onode’s are attached to the physical parent body.

kc::ks_ptr<kc::DataCache<InvDynVectors>> _invdyn_data_cache = nullptr#

Data cache for inverse dynamics data

ATBIDataCaches _atbi_data_caches#

struct with data caches for ATBI data

struct ATBIDataCaches#

ATBI data caches for the body.

Public Members

kc::ks_ptr<kc::DataCache<bool>> matrix_cache = nullptr#

ATBI matrix data cache for the body

kc::ks_ptr<kc::DataCache<bool>> filter_cache = nullptr#

ATBI filter data cache for the body

kc::ks_ptr<kc::DataCache<bool>> smoother_cache = nullptr#

ATBI smoother data cache for the body

kc::ks_ptr<kc::DataCache<bool>> upsilon_matrix_cache = nullptr#

ATBI OSI data cache for the body

struct DumpOptions : public Karana::Core::LockingBase::DumpOptions#

Struct with options for dumpString.

Public Functions

inline DumpOptions &operator=(const DumpOptions &p)#

Assignment operator.

Parameters:

p – The instance to copy from

Returns:

The instance with data copied in.

struct InvDynVectors#

Struct for inverse dynamics interbody spatial force quantities for the embedded bodies.

Public Members

km::Vec f_prime#

the nbodies x 6 raw spatial interaction force from the outboard children bodies together with the local contribution in their pnode frames. We need to process this with Phi_G to get the actual inter-body spatial force for the embedded bodies.

km::Vec f#

the nbodies x 6 interbody spatial force for the embedded bodies. This is Phi_G*f_prime