Description

Class for a generic 3D finite element node, with x,y,z displacement and a direction.

The direction D represents a derivative vector to be used in ANCF elements.

#include <ChNodeFEAxyzD.h>

Inheritance diagram for chrono::fea::ChNodeFEAxyzD:
Collaboration diagram for chrono::fea::ChNodeFEAxyzD:

Public Member Functions

 ChNodeFEAxyzD (ChVector3d initial_pos=VNULL, ChVector3d initial_dir=VECT_X)
 
 ChNodeFEAxyzD (const ChNodeFEAxyzD &other)
 
ChNodeFEAxyzDoperator= (const ChNodeFEAxyzD &other)
 
void SetSlope1 (const ChVector3d &d)
 Set the derivative vector.
 
const ChVector3dGetSlope1 () const
 Get the derivative vector.
 
void SetSlope1Dt (const ChVector3d &dt)
 Set the speed of the derivative vector.
 
const ChVector3dGetSlope1Dt () const
 Get the speed of the derivative vector.
 
void SetSlope1Dt2 (const ChVector3d &dtt)
 Set the acceleration of the derivative vector.
 
const ChVector3dGetSlope1Dt2 () const
 Get the acceleration of the derivative vector.
 
ChVariablesVariablesSlope1 ()
 
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).
 
void SetSlope1Fixed (bool fixed)
 Fix/release the derivative vector states. More...
 
bool IsSlope1Fixed () const
 Return true if the derivative vector states are fixed.
 
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 unsigned int GetNumCoordsPosLevelActive () const override
 Get the actual number of active degrees of freedom.
 
virtual unsigned int GetNumCoordsVelLevelActive () const override
 Get the actual number of active degrees of freedom, derivative.
 
virtual void ArchiveOut (ChArchiveOut &archive) override
 Method to allow serialization of transient data to archives.
 
virtual void ArchiveIn (ChArchiveIn &archive) override
 Method to allow de-serialization of transient data from archives.
 
virtual void SetupInitial (ChSystem *system) override
 Initial setup. Set number of degrees of freedom for this node.
 
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...
 
virtual unsigned int GetLoadableNumCoordsPosLevel () override
 Gets the number of DOFs affected by this element (position part).
 
virtual unsigned int GetLoadableNumCoordsVelLevel () override
 Gets the number of DOFs affected by this element (speed part).
 
virtual void LoadableGetStateBlockPosLevel (int block_offset, ChState &S) override
 Gets all the DOFs packed in a single vector (position part).
 
virtual void LoadableGetStateBlockVelLevel (int block_offset, ChStateDelta &S) override
 Gets all the DOFs packed in a single vector (speed part).
 
virtual void LoadableStateIncrement (const unsigned int off_x, ChState &x_new, const ChState &x, const unsigned int off_v, const ChStateDelta &Dv) override
 Increment all DOFs using a delta.
 
virtual unsigned int GetNumFieldCoords () override
 Number of coordinates in the interpolated field.
 
virtual unsigned int GetSubBlockSize (unsigned int nblock) override
 Get the size of the i-th sub-block of DOFs in global vector.
 
virtual void LoadableGetVariables (std::vector< ChVariables * > &vars) override
 Get the pointers to the contained ChVariables, appending to the mvars vector.
 
virtual void ComputeNF (const double U, const double V, const double W, ChVectorDynamic<> &Qi, double &detJ, const ChVectorDynamic<> &F, ChVectorDynamic<> *state_x, ChVectorDynamic<> *state_w) override
 Evaluate Q = N'*F, for Q generalized lagrangian load, where N is some type of matrix evaluated at point P(U,V,W) assumed in absolute coordinates, and F is a load assumed in absolute coordinates. More...
 
- Public Member Functions inherited from chrono::fea::ChNodeFEAxyz
 ChNodeFEAxyz (ChVector3d initial_pos=VNULL)
 
 ChNodeFEAxyz (const ChNodeFEAxyz &other)
 
ChNodeFEAxyzoperator= (const ChNodeFEAxyz &other)
 
virtual ChVariablesNodeVariables () override
 
virtual void Relax () override
 Set the rest position as the actual position.
 
virtual double GetMass () const override
 Get mass of the node.
 
virtual void SetMass (double m) override
 Set mass of the node.
 
virtual void SetX0 (const ChVector3d &x)
 Set the initial (reference) position.
 
virtual const ChVector3dGetX0 () const
 Get the initial (reference) position.
 
virtual void SetForce (const ChVector3d &frc)
 Set the 3d applied force, in absolute reference.
 
virtual const ChVector3dGetForce () const
 Get the 3d applied force, in absolute reference.
 
virtual ChVariablesGetVariables1 () override
 
- 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.
 
- Public Member Functions inherited from chrono::ChNodeBase
 ChNodeBase (const ChNodeBase &other)
 
ChNodeBaseoperator= (const ChNodeBase &other)
 
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).
 
- Public Member Functions inherited from chrono::ChNodeXYZ
 ChNodeXYZ (const ChVector3d &initial_pos)
 
 ChNodeXYZ (const ChNodeXYZ &other)
 
ChNodeXYZoperator= (const ChNodeXYZ &other)
 
const ChVector3dGetPos () const
 
void SetPos (const ChVector3d &mpos)
 
const ChVector3dGetPosDt () const
 
void SetPosDt (const ChVector3d &mposdt)
 
const ChVector3dGetPosDt2 () const
 
void SetPosDt2 (const ChVector3d &mposdtdt)
 
virtual unsigned int GetNumSubBlocks () override
 Get the number of DOFs sub-blocks.
 
virtual unsigned int GetSubBlockOffset (unsigned int nblock) override
 Get the offset of the specified sub-block of DOFs in global vector.
 
virtual bool IsSubBlockActive (unsigned int nblock) const override
 Check if the specified sub-block of DOFs is active.
 
virtual double GetDensity () override
 This is not needed because not used in quadrature.
 
- Public Member Functions inherited from chrono::ChLoadableUVW
virtual bool IsTetrahedronIntegrationNeeded ()
 If true, use quadrature over u,v,w in [0..1] range as tetrahedron volumetric coords (with z=1-u-v-w) otherwise use default quadrature over u,v,w in [-1..+1] as box isoparametric coords.
 
virtual bool IsTrianglePrismIntegrationNeeded ()
 If true, use quadrature over u,v in [0..1] range as triangle natural coords (with z=1-u-v), and use linear quadrature over w in [-1..+1], otherwise use default quadrature over u,v,w in [-1..+1] as box isoparametric coords.
 

Public Attributes

unsigned int m_dof_actual
 actual number of degrees of freedom
 
- Public Attributes inherited from chrono::fea::ChNodeFEAbase
double m_TotalMass
 Nodal mass obtained from element mass matrix.
 
- Public Attributes inherited from chrono::ChNodeXYZ
ChVector3d pos
 
ChVector3d pos_dt
 
ChVector3d pos_dtdt
 

Protected Member Functions

ChVectorDynamicGetMassDiagonalSlope1 ()
 Get mass of the node (corresponding to the slope derivative).
 

Protected Attributes

ChVariablesGenericDiagonalMassvariables_D
 derivative vector
 
ChVector3d D
 
ChVector3d D_dt
 
ChVector3d D_dtdt
 
- Protected Attributes inherited from chrono::fea::ChNodeFEAxyz
ChVariablesNode variables
 3D node variables, with x,y,z
 
ChVector3d X0
 reference position
 
ChVector3d Force
 applied force
 
- 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 Types inherited from chrono::ChVariableTupleCarrier_1vars< 3 >
typedef ChConstraintTuple_1vars< ChVariableTupleCarrier_1vars< N1 > > type_constraint_tuple
 
- Static Public Attributes inherited from chrono::ChVariableTupleCarrier_1vars< 3 >
static const int nvars1
 

Member Function Documentation

◆ ComputeNF()

void chrono::fea::ChNodeFEAxyzD::ComputeNF ( const double  U,
const double  V,
const double  W,
ChVectorDynamic<> &  Qi,
double &  detJ,
const ChVectorDynamic<> &  F,
ChVectorDynamic<> *  state_x,
ChVectorDynamic<> *  state_w 
)
overridevirtual

Evaluate Q = N'*F, for Q generalized lagrangian load, where N is some type of matrix evaluated at point P(U,V,W) assumed in absolute coordinates, and F is a load assumed in absolute coordinates.

Here, det[J] is unused.

Parameters
Ux coordinate of application point in absolute space
Vy coordinate of application point in absolute space
Wz coordinate of application point in absolute space
QiReturn result of N'*F here, maybe with offset block_offset
detJReturn det[J] here
FInput F vector, containing Force xyz in absolute coords and a 'pseudo' torque.
state_xif != 0, update state (pos. part) to this, then evaluate Q
state_wif != 0, update state (speed part) to this, then evaluate Q

Reimplemented from chrono::ChNodeXYZ.

Reimplemented in chrono::fea::ChNodeFEAxyzDDD, and chrono::fea::ChNodeFEAxyzDD.

◆ SetFixed()

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

Fix/release this node.

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

Reimplemented from chrono::fea::ChNodeFEAxyz.

Reimplemented in chrono::fea::ChNodeFEAxyzDDD, and chrono::fea::ChNodeFEAxyzDD.

◆ SetSlope1Fixed()

void chrono::fea::ChNodeFEAxyzD::SetSlope1Fixed ( bool  fixed)

Fix/release the derivative vector states.

If fixed, these states are not changed by the solver.

◆ VariablesFbIncrementMq()

void chrono::fea::ChNodeFEAxyzD::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::fea::ChNodeFEAxyz.

Reimplemented in chrono::fea::ChNodeFEAxyzDDD, and chrono::fea::ChNodeFEAxyzDD.

◆ VariablesFbLoadForces()

void chrono::fea::ChNodeFEAxyzD::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::fea::ChNodeFEAxyz.

Reimplemented in chrono::fea::ChNodeFEAxyzDDD, and chrono::fea::ChNodeFEAxyzDD.

◆ VariablesQbIncrementPosition()

void chrono::fea::ChNodeFEAxyzD::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::fea::ChNodeFEAxyz.

Reimplemented in chrono::fea::ChNodeFEAxyzDDD, and chrono::fea::ChNodeFEAxyzDD.

◆ VariablesQbSetSpeed()

void chrono::fea::ChNodeFEAxyzD::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::fea::ChNodeFEAxyz.

Reimplemented in chrono::fea::ChNodeFEAxyzDDD, and chrono::fea::ChNodeFEAxyzDD.


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