chrono::ChLoadXYZnodeXYZnode Class Referenceabstract

Description

Base class for loads representing a concentrated force acting between two ChNodeXYZ.

Users should inherit from this and implement a custom ComputeForce(), this is enough to have the load working. Note: there are already some predefined examples, ie. children classes with common concrete implementations, such as ChLoadXYZnodeXYZnodeSpring ChLoadXYZnodeXYZnodeBushing etc.; take inspiration from those example to inherit your own load if they are not enough.

#include <ChLoadsXYZnode.h>

Inheritance diagram for chrono::ChLoadXYZnodeXYZnode:
Collaboration diagram for chrono::ChLoadXYZnodeXYZnode:

Public Member Functions

 ChLoadXYZnodeXYZnode (std::shared_ptr< ChNodeXYZ > mnodeA, std::shared_ptr< ChNodeXYZ > mnodeB)
 
virtual void ComputeForce (const ChVector<> &rel_pos, const ChVector<> &rel_vel, ChVector<> &abs_force)=0
 Compute the force on the node A, in absolute coordsystem, given position and speed of node respect to the other. More...
 
virtual void ComputeQ (ChState *state_x, ChStateDelta *state_w) override
 Compute Q, the generalized load. More...
 
ChVector GetForce () const
 For diagnosis purposes, this can return the actual last computed value of the applied force, expressed in absolute coordinate system, assumed applied to node.
 
- Public Member Functions inherited from chrono::ChLoadCustomMultiple
 ChLoadCustomMultiple (std::vector< std::shared_ptr< ChLoadable >> &mloadables)
 
 ChLoadCustomMultiple (std::shared_ptr< ChLoadable > mloadableA, std::shared_ptr< ChLoadable > mloadableB)
 
 ChLoadCustomMultiple (std::shared_ptr< ChLoadable > mloadableA, std::shared_ptr< ChLoadable > mloadableB, std::shared_ptr< ChLoadable > mloadableC)
 
virtual int LoadGet_ndof_x () override
 Gets the number of DOFs affected by this load (position part)
 
virtual int LoadGet_ndof_w () override
 Gets the number of DOFs affected by this load (speed part)
 
virtual void LoadGetStateBlock_x (ChState &mD) override
 Gets all the current DOFs packed in a single vector (position part)
 
virtual void LoadGetStateBlock_w (ChStateDelta &mD) override
 Gets all the current DOFs packed in a single vector (speed part)
 
virtual void LoadStateIncrement (const ChState &x, const ChStateDelta &dw, ChState &x_new) override
 Increment a packed state (ex. More...
 
virtual int LoadGet_field_ncoords () override
 Number of coordinates in the interpolated field, ex=3 for a tetrahedron finite element or a cable, = 1 for a thermal problem, etc.
 
virtual void ComputeJacobian (ChState *state_x, ChStateDelta *state_w, ChMatrixRef mK, ChMatrixRef mR, ChMatrixRef mM) override
 Compute jacobians (default fallback). More...
 
virtual void LoadIntLoadResidual_F (ChVectorDynamic<> &R, double c) override
 Adds the internal loads Q (pasted at global nodes offsets) into a global vector R, multiplied by a scaling factor c, as R += forces * c.
 
virtual void LoadIntLoadResidual_Mv (ChVectorDynamic<> &R, const ChVectorDynamic<> &w, double c) override
 Default fallback: compute jacobians via ComputeJacobian(), then use M=-dQ/da to do R += c*M*w. More...
 
virtual void CreateJacobianMatrices () override
 Create the jacobian loads if needed, and also set the ChVariables referenced by the sparse KRM block.
 
virtual ChVectorDynamicGetQ ()
 Access the generalized load vector Q.
 
- Public Member Functions inherited from chrono::ChLoadBase
ChLoadJacobiansGetJacobians ()
 Access the jacobians (if any, i.e. if this is a stiff load)
 
virtual void InjectKRMmatrices (ChSystemDescriptor &mdescriptor)
 Tell to a system descriptor that there are item(s) of type ChKblock in this object (for further passing it to a solver)
 
virtual void KRMmatricesLoad (double Kfactor, double Rfactor, double Mfactor)
 Adds the current stiffness K and damping R and mass M matrices in encapsulated ChKblock item(s), if any. More...
 
- Public Member Functions inherited from chrono::ChObj
 ChObj (const ChObj &other)
 
virtual ChObjClone () const =0
 "Virtual" copy constructor. More...
 
int GetIdentifier () const
 Gets the numerical identifier of the object.
 
void SetIdentifier (int id)
 Sets the numerical identifier of the object.
 
double GetChTime () const
 Gets the simulation time of this object.
 
void SetChTime (double m_time)
 Sets the simulation time of this object.
 
const char * GetName () const
 Gets the name of the object as C Ascii null-terminated string -for reading only!
 
void SetName (const char myname[])
 Sets the name of this object, as ascii string.
 
std::string GetNameString () const
 Gets the name of the object as C Ascii null-terminated string.
 
void SetNameString (const std::string &myname)
 Sets the name of this object, as std::string.
 
void MFlagsSetAllOFF (int &mflag)
 
void MFlagsSetAllON (int &mflag)
 
void MFlagSetON (int &mflag, int mask)
 
void MFlagSetOFF (int &mflag, int mask)
 
int MFlagGet (int &mflag, int mask)
 
virtual void ArchiveOUT (ChArchiveOut &marchive)
 Method to allow serialization of transient data to archives.
 
virtual void ArchiveIN (ChArchiveIn &marchive)
 Method to allow de-serialization of transient data from archives.
 
virtual std::string & ArchiveContainerName ()
 

Protected Member Functions

virtual bool IsStiff () override
 Inherited classes could override this and return true, if the load benefits from a jacobian when using implicit integrators. More...
 
virtual void Update (double time) override
 Update: this is called at least at each time step. More...
 

Protected Attributes

ChVector computed_abs_force
 
- Protected Attributes inherited from chrono::ChLoadBase
ChLoadJacobiansjacobians
 
- Protected Attributes inherited from chrono::ChObj
double ChTime
 the time of simulation for the object
 

Additional Inherited Members

- Public Attributes inherited from chrono::ChLoadCustomMultiple
std::vector< std::shared_ptr< ChLoadable > > loadables
 
ChVectorDynamic load_Q
 

Constructor & Destructor Documentation

◆ ChLoadXYZnodeXYZnode()

chrono::ChLoadXYZnodeXYZnode::ChLoadXYZnodeXYZnode ( std::shared_ptr< ChNodeXYZ mnodeA,
std::shared_ptr< ChNodeXYZ mnodeB 
)
Parameters
mnodeAnode A to apply load to
mnodeBnode B to apply load to, as reaction

Member Function Documentation

◆ ComputeForce()

virtual void chrono::ChLoadXYZnodeXYZnode::ComputeForce ( const ChVector<> &  rel_pos,
const ChVector<> &  rel_vel,
ChVector<> &  abs_force 
)
pure virtual

Compute the force on the node A, in absolute coordsystem, given position and speed of node respect to the other.

Node B gets the opposite. Inherited classes MUST IMPLEMENT THIS.

Implemented in chrono::ChLoadXYZnodeXYZnodeBushing, and chrono::ChLoadXYZnodeXYZnodeSpring.

◆ ComputeQ()

void chrono::ChLoadXYZnodeXYZnode::ComputeQ ( ChState state_x,
ChStateDelta state_w 
)
overridevirtual

Compute Q, the generalized load.

Called automatically at each Update().

Parameters
state_xstate position to evaluate Q
state_wstate speed to evaluate Q

Implements chrono::ChLoadBase.

◆ IsStiff()

virtual bool chrono::ChLoadXYZnodeXYZnode::IsStiff ( )
inlineoverrideprotectedvirtual

Inherited classes could override this and return true, if the load benefits from a jacobian when using implicit integrators.


Implements chrono::ChLoadBase.

Reimplemented in chrono::ChLoadXYZnodeXYZnodeBushing, and chrono::ChLoadXYZnodeXYZnodeSpring.

◆ Update()

void chrono::ChLoadXYZnodeXYZnode::Update ( double  time)
overrideprotectedvirtual

Update: this is called at least at each time step.

  • It recomputes the generalized load Q vector(s)
  • It recomputes the jacobian(s) K,R,M in case of stiff load Q and jacobians assumed evaluated at the current state. Jacobian structures are automatically allocated if needed.

Reimplemented from chrono::ChLoadBase.


The documentation for this class was generated from the following files:
  • /builds/uwsbel/chrono/src/chrono/physics/ChLoadsXYZnode.h
  • /builds/uwsbel/chrono/src/chrono/physics/ChLoadsXYZnode.cpp