Description

Generic finite element node with 9 degrees of freedom representing curvature.

#include <ChNodeFEAcurv.h>

Inheritance diagram for chrono::fea::ChNodeFEAcurv:
Collaboration diagram for chrono::fea::ChNodeFEAcurv:

Public Member Functions

 ChNodeFEAcurv (const ChVector3d &rxx=VNULL, const ChVector3d &ryy=VNULL, const ChVector3d &rzz=VNULL)
 
 ChNodeFEAcurv (const ChNodeFEAcurv &other)
 
ChNodeFEAcurvoperator= (const ChNodeFEAcurv &other)
 
void SetCurvatureXX (const ChVector3d &rxx)
 Set the xx 2nd derivative of position vector.
 
const ChVector3dGetCurvatureXX () const
 Get the xx 2nd derivative of position vector.
 
void SetCurvatureYY (const ChVector3d &ryy)
 Set the yy 2nd derivative of position vector.
 
const ChVector3dGetCurvatureYY () const
 Get the yy 2nd derivative of position vector.
 
void SetCurvatureZZ (const ChVector3d &rzz)
 Set the zz 2nd derivative of position vector.
 
const ChVector3dGetCurvatureZZ () const
 Get the zz 2nd derivative of position vector.
 
void SetCurvatureXX_dt (const ChVector3d &rxx_dt)
 Set the time derivative of the xx 2nd derivative of position vector.
 
const ChVector3dGetCurvatureXX_dt () const
 Get the time derivative of the xx 2nd derivative of position vector.
 
void SetCurvatureYY_dt (const ChVector3d &ryy_dt)
 Set the time derivative of the yy 2nd derivative of position vector.
 
const ChVector3dGetCurvatureYY_dt () const
 Get the time derivative of the yy 2nd derivative of position vector.
 
void SetCurvatureZZ_dt (const ChVector3d &rzz_dt)
 Set the time derivative of the zz 2nd derivative of position vector.
 
const ChVector3dGetCurvatureZZ_dt () const
 Get the time derivative of the zz 2nd derivative of position vector.
 
ChVectorDynamicGetMassDiagonal ()
 Get mass of the node. / TODO is this even meaningful/needed for this type of node?
 
void SetMass (double mass)
 Set mass of the node. / TODO is this even meaningful/needed for this type of node?
 
ChVariablesVariables ()
 
virtual void Relax () override
 Reset the 2nd derivatives of position vector and their time derivatives.
 
virtual void ForceToRest () override
 Reset to no speed and acceleration.
 
virtual void SetFixed (bool fixed) override
 Fix/release this node. More...
 
virtual bool IsFixed () const override
 Return true if the node is fixed (i.e., its state variables are not changed by the solver).
 
virtual unsigned int GetNumCoordsPosLevel () const override
 Get the number of degrees of freedom.
 
virtual unsigned int GetNumCoordsVelLevel () const override
 Get the number of degrees of freedom, derivative.
 
virtual void ArchiveOut (ChArchiveOut &archive_out) override
 Method to allow serialization of transient data to archives.
 
virtual void ArchiveIn (ChArchiveIn &archive_in) override
 Method to allow de-serialization of transient data from archives.
 
virtual void NodeIntStateGather (const unsigned int off_x, ChState &x, const unsigned int off_v, ChStateDelta &v, double &T) override
 
virtual void NodeIntStateScatter (const unsigned int off_x, const ChState &x, const unsigned int off_v, const ChStateDelta &v, const double T) override
 
virtual void NodeIntStateGatherAcceleration (const unsigned int off_a, ChStateDelta &a) override
 
virtual void NodeIntStateScatterAcceleration (const unsigned int off_a, const ChStateDelta &a) override
 
virtual void NodeIntStateIncrement (const unsigned int off_x, ChState &x_new, const ChState &x, const unsigned int off_v, const ChStateDelta &Dv) override
 
virtual void NodeIntStateGetIncrement (const unsigned int off_x, const ChState &x_new, const ChState &x, const unsigned int off_v, ChStateDelta &Dv) override
 
virtual void NodeIntLoadResidual_F (const unsigned int off, ChVectorDynamic<> &R, const double c) override
 
virtual void NodeIntLoadResidual_Mv (const unsigned int off, ChVectorDynamic<> &R, const ChVectorDynamic<> &w, const double c) override
 
virtual void NodeIntLoadLumpedMass_Md (const unsigned int off, ChVectorDynamic<> &Md, double &error, const double c) override
 
virtual void NodeIntToDescriptor (const unsigned int off_v, const ChStateDelta &v, const ChVectorDynamic<> &R) override
 
virtual void NodeIntFromDescriptor (const unsigned int off_v, ChStateDelta &v) override
 
virtual void InjectVariables (ChSystemDescriptor &descriptor) override
 Register with the given system descriptor any ChVariable objects associated with this item.
 
virtual void VariablesFbReset () override
 Set the 'fb' part (the known term) of the encapsulated ChVariables to zero.
 
virtual void VariablesFbLoadForces (double factor=1) override
 Add the current forces (applied to node) into the encapsulated ChVariables. More...
 
virtual void VariablesQbLoadSpeed () override
 Initialize the 'qb' part of the ChVariables with the current value of speeds.
 
virtual void VariablesQbSetSpeed (double step=0) override
 Fetch the item speed (ex. More...
 
virtual void VariablesFbIncrementMq () override
 Add M*q (masses multiplied current 'qb') to Fb, ex. More...
 
virtual void VariablesQbIncrementPosition (double step) override
 Increment node positions by the 'qb' part of the ChVariables, multiplied by a 'step' factor. More...
 
- Public Member Functions inherited from chrono::fea::ChNodeFEAbase
virtual void SetIndex (unsigned int mindex)
 Sets the global index of the node.
 
virtual unsigned int GetIndex ()
 Gets the global index of the node.
 
virtual void SetupInitial (ChSystem *system)
 Initial setup.
 
- Public Member Functions inherited from chrono::ChNodeBase
 ChNodeBase (const ChNodeBase &other)
 
ChNodeBaseoperator= (const ChNodeBase &other)
 
virtual unsigned int GetNumCoordsPosLevelActive () const
 Get the actual number of active degrees of freedom. More...
 
virtual unsigned int GetNumCoordsVelLevelActive () const
 Get the actual number of active degrees of freedom, derivative. More...
 
virtual bool IsAllCoordsActive () const
 Return true if all node DOFs are active (no node variable is fixed).
 
unsigned int NodeGetOffsetPosLevel ()
 Get offset in the state vector (position part).
 
unsigned int NodeGetOffsetVelLevel ()
 Get offset in the state vector (speed part).
 
void NodeSetOffsetPosLevel (const unsigned int moff)
 Set offset in the state vector (position part).
 
void NodeSetOffsetVelLevel (const unsigned int moff)
 Set offset in the state vector (speed part).
 

Protected Attributes

ChVariablesGenericDiagonalMassm_variables
 
ChVector3d m_rxx
 
ChVector3d m_ryy
 
ChVector3d m_rzz
 
ChVector3d m_rxx_dt
 
ChVector3d m_ryy_dt
 
ChVector3d m_rzz_dt
 
ChVector3d m_rxx_dtdt
 
ChVector3d m_ryy_dtdt
 
ChVector3d m_rzz_dtdt
 
- Protected Attributes inherited from chrono::fea::ChNodeFEAbase
unsigned int g_index
 global node index
 
- Protected Attributes inherited from chrono::ChNodeBase
unsigned int offset_x
 offset in vector of state (position part)
 
unsigned int offset_w
 offset in vector of state (speed part)
 

Additional Inherited Members

- Public Attributes inherited from chrono::fea::ChNodeFEAbase
double m_TotalMass
 Nodal mass obtained from element mass matrix.
 

Constructor & Destructor Documentation

◆ ChNodeFEAcurv()

chrono::fea::ChNodeFEAcurv::ChNodeFEAcurv ( const ChVector3d rxx = VNULL,
const ChVector3d ryy = VNULL,
const ChVector3d rzz = VNULL 
)
Parameters
rxxinitial value of xx 2nd derivative of position vector
ryyinitial value of yy 2nd derivative of position vector
rzzinitial value of zz 2nd derivative of position vector

Member Function Documentation

◆ SetFixed()

void chrono::fea::ChNodeFEAcurv::SetFixed ( bool  fixed)
overridevirtual

Fix/release this node.

If fixed, its state variables are not changed by the solver.

Implements chrono::fea::ChNodeFEAbase.

◆ VariablesFbIncrementMq()

void chrono::fea::ChNodeFEAcurv::VariablesFbIncrementMq ( )
overridevirtual

Add M*q (masses multiplied current 'qb') to Fb, ex.

if qb is initialized with v_old using VariablesQbLoadSpeed, this method can be used in timestepping schemes that do: M*v_new = M*v_old + forces*dt

Reimplemented from chrono::ChNodeBase.

◆ VariablesFbLoadForces()

void chrono::fea::ChNodeFEAcurv::VariablesFbLoadForces ( double  factor = 1)
overridevirtual

Add the current forces (applied to node) into the encapsulated ChVariables.

Include in the 'fb' part: qf+=forces*factor

Reimplemented from chrono::ChNodeBase.

◆ VariablesQbIncrementPosition()

void chrono::fea::ChNodeFEAcurv::VariablesQbIncrementPosition ( double  step)
overridevirtual

Increment node positions by the 'qb' part of the ChVariables, multiplied by a 'step' factor.

pos+=qb*step If qb is a speed, this behaves like a single step of 1-st order numerical integration (Euler integration).

Reimplemented from chrono::ChNodeBase.

◆ VariablesQbSetSpeed()

void chrono::fea::ChNodeFEAcurv::VariablesQbSetSpeed ( double  step = 0)
overridevirtual

Fetch the item speed (ex.

linear velocity, in xyz nodes) from the 'qb' part of the ChVariables and sets it as the current item speed. If 'step' is not 0, also should compute the approximate acceleration of the item using backward differences, that is accel=(new_speed-old_speed)/step. Mostly used after the solver provided the solution in ChVariables.

Reimplemented from chrono::ChNodeBase.


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