Description
Class for a generic 3D finite element node, with x,y,z displacement and a direction.
The direction D represents a gradient vector to be used in ANCF elements.
#include <ChNodeFEAxyzD.h>
Public Member Functions | |
ChNodeFEAxyzD (ChVector<> initial_pos=VNULL, ChVector<> initial_dir=VECT_X) | |
ChNodeFEAxyzD (const ChNodeFEAxyzD &other) | |
ChNodeFEAxyzD & | operator= (const ChNodeFEAxyzD &other) |
void | SetD (ChVector<> mD) |
Set the direction. | |
const ChVector & | GetD () const |
Get the direction. | |
void | SetD_dt (ChVector<> mD) |
Set the direction speed. | |
const ChVector & | GetD_dt () const |
Get the direction speed. | |
void | SetD_dtdt (ChVector<> mD) |
Set the direction acceleration. | |
const ChVector & | GetD_dtdt () const |
Get the direction acceleration. | |
ChVariables & | Variables_D () |
virtual void | SetNoSpeedNoAcceleration () override |
Reset to no speed and acceleration. | |
virtual ChVectorDynamic & | GetMassDiagonal () |
Get mass of the node. | |
virtual void | SetFixed (bool mev) override |
Sets the 'fixed' state of the node. More... | |
virtual bool | GetFixed () override |
Gets the 'fixed' state of the node. | |
virtual int | Get_ndof_x () const override |
Get the number of degrees of freedom. | |
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 | 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 | 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 &mdescriptor) override |
Tell to a system descriptor that there are variables of type ChVariables in this object (for further passing it to a solver) | |
virtual void | VariablesFbReset () override |
Sets the 'fb' part (the known term) of the encapsulated ChVariables to zero. | |
virtual void | VariablesFbLoadForces (double factor=1) override |
Adds the current forces (applied to node) into the encapsulated ChVariables, in the 'fb' part: qf+=forces*factor. | |
virtual void | VariablesQbLoadSpeed () override |
Initialize the 'qb' part of the ChVariables with the current value of speeds. | |
virtual void | VariablesQbSetSpeed (double step=0) override |
Fetches the item speed (ex. More... | |
virtual void | VariablesFbIncrementMq () override |
Adds 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 int | LoadableGet_ndof_x () override |
Gets the number of DOFs affected by this element (position part) | |
virtual int | LoadableGet_ndof_w () override |
Gets the number of DOFs affected by this element (speed part) | |
virtual void | LoadableGetStateBlock_x (int block_offset, ChState &mD) override |
Gets all the DOFs packed in a single vector (position part) | |
virtual void | LoadableGetStateBlock_w (int block_offset, ChStateDelta &mD) 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 int | Get_field_ncoords () override |
Number of coordinates in the interpolated field, ex=3 for a tetrahedron finite element or a cable, etc. More... | |
virtual unsigned int | GetSubBlockSize (int nblock) override |
Get the size of the i-th sub-block of DOFs in global vector. | |
virtual void | LoadableGetVariables (std::vector< ChVariables * > &mvars) 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... | |
virtual void | ArchiveOUT (ChArchiveOut &marchive) override |
Method to allow serialization of transient data to archives. | |
virtual void | ArchiveIN (ChArchiveIn &marchive) override |
Method to allow de serialization of transient data from archives. | |
Public Member Functions inherited from chrono::fea::ChNodeFEAxyz | |
ChNodeFEAxyz (ChVector<> initial_pos=VNULL) | |
ChNodeFEAxyz (const ChNodeFEAxyz &other) | |
ChNodeFEAxyz & | operator= (const ChNodeFEAxyz &other) |
virtual ChVariablesNode & | Variables () 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 mm) override |
Set mass of the node. | |
virtual void | SetX0 (ChVector<> mx) |
Set the initial (reference) position. | |
virtual ChVector | GetX0 () |
Get the initial (reference) position. | |
virtual void | SetForce (ChVector<> mf) |
Set the 3d applied force, in absolute reference. | |
virtual ChVector | GetForce () |
Get the 3d applied force, in absolute reference. | |
virtual ChVariables * | GetVariables1 () 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) | |
ChNodeBase & | operator= (const ChNodeBase &other) |
virtual int | Get_ndof_w () const |
Get the number of degrees of freedom, derivative This might be different from ndof if quaternions are used for rotations, as derivative might be angular velocity. | |
unsigned int | NodeGetOffset_x () |
Get offset in the state vector (position part) | |
unsigned int | NodeGetOffset_w () |
Get offset in the state vector (speed part) | |
void | NodeSetOffset_x (const unsigned int moff) |
Set offset in the state vector (position part) | |
void | NodeSetOffset_w (const unsigned int moff) |
Set offset in the state vector (speed part) | |
Public Member Functions inherited from chrono::ChNodeXYZ | |
ChNodeXYZ (const ChVector<> &initial_pos) | |
ChNodeXYZ (const ChNodeXYZ &other) | |
ChNodeXYZ & | operator= (const ChNodeXYZ &other) |
const ChVector & | GetPos () const |
void | SetPos (const ChVector<> &mpos) |
const ChVector & | GetPos_dt () const |
void | SetPos_dt (const ChVector<> &mposdt) |
const ChVector & | GetPos_dtdt () const |
void | SetPos_dtdt (const ChVector<> &mposdtdt) |
virtual int | GetSubBlocks () override |
Tell the number of DOFs blocks (ex. =1 for a body, =4 for a tetrahedron, etc.) | |
virtual unsigned int | GetSubBlockOffset (int nblock) override |
Get the offset of the i-th sub-block of DOFs in global vector. | |
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 | |
ChVector | D |
ChVector | D_dt |
ChVector | D_dtdt |
Public Attributes inherited from chrono::fea::ChNodeFEAbase | |
double | m_TotalMass |
Nodal mass obtained from element masss matrix. | |
Public Attributes inherited from chrono::ChNodeXYZ | |
ChVector | pos |
ChVector | pos_dt |
ChVector | pos_dtdt |
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 |
Protected Attributes inherited from chrono::fea::ChNodeFEAxyz | |
ChVariablesNode | variables |
ChVector | X0 |
3D node variables, with x,y,z More... | |
ChVector | 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) | |
Member Function Documentation
◆ ComputeNF()
|
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.
The det[J] is unused.
- Parameters
-
U x coordinate of application point in absolute space V y coordinate of application point in absolute space W z coordinate of application point in absolute space Qi Return result of N'*F here, maybe with offset block_offset detJ Return det[J] here F Input F vector, containing Force xyz in absolute coords and a 'pseudo' torque. state_x if != 0, update state (pos. part) to this, then evaluate Q state_w if != 0, update state (speed part) to this, then evaluate Q
Reimplemented from chrono::ChNodeXYZ.
Reimplemented in chrono::fea::ChNodeFEAxyzDD.
◆ Get_field_ncoords()
|
inlineoverridevirtual |
Number of coordinates in the interpolated field, ex=3 for a tetrahedron finite element or a cable, etc.
Here is 6: xyz displ + xyz rots
Reimplemented from chrono::ChNodeXYZ.
Reimplemented in chrono::fea::ChNodeFEAxyzDD.
◆ SetFixed()
|
overridevirtual |
Sets the 'fixed' state of the node.
If true, it does not move respect to the absolute world, despite constraints, forces, etc.
Reimplemented from chrono::fea::ChNodeFEAxyz.
Reimplemented in chrono::fea::ChNodeFEAxyzDD.
◆ VariablesFbIncrementMq()
|
overridevirtual |
Adds 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::ChNodeFEAxyzDD.
◆ VariablesQbIncrementPosition()
|
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 (Eulero integration).
Reimplemented from chrono::fea::ChNodeFEAxyz.
Reimplemented in chrono::fea::ChNodeFEAxyzDD.
◆ VariablesQbSetSpeed()
|
overridevirtual |
Fetches 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::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