Description
Class for a generic ED finite element node, with x,y,z displacement and a 3D rotation.
This is the typical node that can be used for beams, etc.
#include <ChNodeFEAxyzrot.h>
Public Member Functions | |
ChNodeFEAxyzrot (ChFrame<> initialf=ChFrame<>()) | |
ChNodeFEAxyzrot (const ChNodeFEAxyzrot &other) | |
ChNodeFEAxyzrot & | operator= (const ChNodeFEAxyzrot &other) |
virtual ChVariables & | Variables () override |
Return a reference to the encapsulated variables, representing states (pos, speed or accel.) and forces. | |
virtual void | Relax () override |
Set the rest position as the actual position. | |
virtual void | SetNoSpeedNoAcceleration () override |
Reset to no speed and acceleration. | |
virtual void | SetFixed (bool mev) override |
Set the 'fixed' state of the node. More... | |
virtual bool | GetFixed () override |
Get the 'fixed' state of the node. More... | |
double | GetMass () |
Get atomic mass of the node. | |
void | SetMass (double mm) |
Set atomic mass of the node. | |
ChMatrix33 & | GetInertia () |
Access atomic inertia of the node. | |
void | SetX0 (ChFrame<> mx) |
Set the initial (reference) frame. | |
const ChFrame & | GetX0 () const |
Get the initial (reference) frame. | |
ChFrame & | GetX0ref () |
Access the initial (reference) frame. | |
void | SetForce (ChVector<> mf) |
Set the 3d applied force, in absolute reference. | |
const ChVector & | GetForce () const |
Get the 3d applied force, in absolute reference. | |
void | SetTorque (ChVector<> mf) |
Set the 3d applied torque, in node reference. | |
const ChVector & | GetTorque () const |
Get the 3d applied torque, in node reference. | |
ChFrameMoving & | Frame () |
Access the frame of the node - in absolute csys, with infos on actual position, speed, acceleration, etc. | |
virtual int | Get_ndof_x () const override |
Get the number of degrees of freedom (7 because quaternion for rotation). | |
virtual int | Get_ndof_w () const override |
Get the number of degrees of freedom, derivative (6 because angular velocity for rotation derivative). | |
virtual ChVariables * | GetVariables1 () override |
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 int | GetSubBlocks () override |
Get the number of DOFs sub-blocks. | |
virtual unsigned int | GetSubBlockOffset (int nblock) override |
Get the offset of the specified sub-block of DOFs in global vector. | |
virtual unsigned int | GetSubBlockSize (int nblock) override |
Get the size of the specified sub-block of DOFs in global vector. | |
virtual bool | IsSubBlockActive (int nblock) const override |
Check if the specified sub-block of DOFs is active. | |
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 double | GetDensity () override |
This is not needed because not used in quadrature. | |
virtual void | ArchiveOUT (ChArchiveOut &marchive) override |
Method to allow serialization of transient data to archives. | |
virtual void | ArchiveIN (ChArchiveIn &marchive) override |
Method to allow deserialization of transient data from archives. | |
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) |
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::ChBodyFrame | |
ChBodyFrame (const ChBodyFrame &other) | |
void | To_abs_forcetorque (const ChVector<> &force, const ChVector<> &appl_point, bool local, ChVector<> &resultforce, ChVector<> &resulttorque) |
Transform a force applied at a point on the body into a force and moment applied to the COM and expressed in the absolute frame. More... | |
Public Member Functions inherited from chrono::ChFrameMoving< double > | |
ChFrameMoving (const ChVector< double > &mv=ChVector< double >(0, 0, 0), const ChQuaternion< double > &mq=ChQuaternion< double >(1, 0, 0, 0)) | |
Construct from pos and rot (as a quaternion) | |
ChFrameMoving (const ChVector< double > &mv, const ChMatrix33< double > &ma) | |
Construct from pos and rotation (as a 3x3 matrix) | |
ChFrameMoving (const ChCoordsys< double > &mc) | |
Construct from a coordsys. | |
ChFrameMoving (const ChFrame< double > &mc) | |
Construct from a frame. | |
ChFrameMoving (const ChFrameMoving< double > &other) | |
Copy constructor, build from another moving frame. | |
virtual | ~ChFrameMoving () |
Destructor. | |
ChFrameMoving< double > & | operator= (const ChFrameMoving< double > &other) |
Assignment operator: copy from another moving frame. | |
ChFrameMoving< double > & | operator= (const ChFrame< double > &other) |
Assignment operator: copy from another frame. | |
bool | operator== (const ChFrameMoving< double > &other) const |
Returns true for identical frames. | |
bool | operator!= (const ChFrameMoving< double > &other) const |
Returns true for different frames. | |
ChFrameMoving< double > | operator>> (const ChFrameMoving< double > &Fb) const |
The '>>' operator transforms a coordinate system, so transformations can be represented with this syntax: new_frame = old_frame >> tr_frame; For a sequence of transformations, i.e. More... | |
ChFrameMoving< double > | operator* (const ChFrameMoving< double > &Fb) const |
The '*' operator transforms a coordinate system, so transformations can be represented with this syntax: new_frame = tr_frame * old_frame; For a sequence of transformations, i.e. More... | |
ChFrameMoving< double > & | operator>>= (const ChFrameMoving< double > &T) |
Performs pre-multiplication of this frame by another frame, for example: A>>=T means A'=T*A ; or A'=A >> T. | |
ChFrameMoving< double > & | operator>>= (const ChVector< double > &D) |
Performs pre-multiplication of this frame by a vector D, to 'move' by a displacement D: | |
ChFrameMoving< double > & | operator>>= (const ChQuaternion< double > &R) |
Performs pre-multiplication of this frame by a quaternion R, to 'rotate' it by R: | |
ChFrameMoving< double > & | operator>>= (const ChCoordsys< double > &F) |
Performs pre-multiplication of this frame by a ChCoordsys F: | |
ChFrameMoving< double > & | operator>>= (const ChFrame< double > &F) |
Performs pre-multiplication of this frame by a ChFrame F: | |
ChFrameMoving< double > & | operator%= (const ChFrameMoving< double > &T) |
Performs pre-multiplication of this frame by another frame, for example: A%=T means A'=T*A ; or A'=A >> T Note: DEPRECATED, use >>= instead. | |
ChFrameMoving< double > & | operator*= (const ChFrameMoving< double > &T) |
Performs post-multiplication of this frame by another frame, for example: A*=T means A'=A*T ; or A'=T >> A. | |
ChCoordsys< double > & | GetCoord_dt () |
Return both current rotation and translation speeds as a coordsystem object, with vector and quaternion. | |
const ChCoordsys< double > & | GetCoord_dt () const |
ChCoordsys< double > & | GetCoord_dtdt () |
Return both current rotation and translation accelerations as a coordsystem object, with vector and quaternion. | |
const ChCoordsys< double > & | GetCoord_dtdt () const |
ChVector< double > & | GetPos_dt () |
Return the current speed as a 3d vector. | |
const ChVector< double > & | GetPos_dt () const |
ChVector< double > & | GetPos_dtdt () |
Return the current acceleration as a 3d vector. | |
const ChVector< double > & | GetPos_dtdt () const |
ChQuaternion< double > & | GetRot_dt () |
Return the current rotation speed as a quaternion. | |
const ChQuaternion< double > & | GetRot_dt () const |
ChQuaternion< double > & | GetRot_dtdt () |
Return the current rotation acceleration as a quaternion. | |
const ChQuaternion< double > & | GetRot_dtdt () const |
ChVector< double > | GetWvel_loc () const |
Computes the actual angular speed (expressed in local coords) | |
ChVector< double > | GetWvel_par () const |
Computes the actual angular speed (expressed in parent coords) | |
ChVector< double > | GetWacc_loc () const |
Computes the actual angular acceleration (expressed in local coords) | |
ChVector< double > | GetWacc_par () const |
Computes the actual angular acceleration (expressed in parent coords) | |
virtual void | SetCoord_dt (const ChCoordsys< double > &mcoord_dt) |
Set both linear speed and rotation speed as a single ChCoordsys derivative. | |
virtual void | SetPos_dt (const ChVector< double > &mvel) |
Set the linear speed. | |
virtual void | SetRot_dt (const ChQuaternion< double > &mrot_dt) |
Set the rotation speed as a quaternion. More... | |
virtual void | SetWvel_loc (const ChVector< double > &wl) |
Set the rotation speed from given angular speed (expressed in local csys) | |
virtual void | SetWvel_par (const ChVector< double > &wp) |
Set the rotation speed from given angular speed (expressed in parent csys) | |
virtual void | SetCoord_dtdt (const ChCoordsys< double > &mcoord_dtdt) |
Set both linear acceleration and rotation acceleration as a single ChCoordsys derivative. | |
virtual void | SetPos_dtdt (const ChVector< double > &macc) |
Set the linear acceleration. | |
virtual void | SetRot_dtdt (const ChQuaternion< double > &mrot_dtdt) |
Set the rotation acceleration as a quaternion derivative. More... | |
virtual void | SetWacc_loc (const ChVector< double > &al) |
Set the rotation acceleration from given angular acceleration (expressed in local csys) Note: even when the local angular acceleration is zero, you are still encouraged to call this method bacause q_dtdt might be nonzero due to nonzero q_dt in case of rotational motion. | |
virtual void | SetWacc_par (const ChVector< double > &ap) |
Set the rotation speed from given angular speed (expressed in parent csys) | |
void | Compute_Adt (ChMatrix33< double > &mA_dt) const |
Computes the time derivative of rotation matrix, mAdt. | |
void | Compute_Adtdt (ChMatrix33< double > &mA_dtdt) |
Computes the 2nd time derivative of rotation matrix, mAdtdt. | |
ChMatrix33< double > | GetA_dt () |
Computes and returns an Adt matrix (-note: prefer using Compute_Adt() directly for better performance) | |
ChMatrix33< double > | GetA_dtdt () |
Computes and returns an Adt matrix (-note: prefer using Compute_Adtdt() directly for better performance) | |
void | ConcatenatePreTransformation (const ChFrameMoving< double > &T) |
Apply a transformation (rotation and translation) represented by another ChFrameMoving T. More... | |
void | ConcatenatePostTransformation (const ChFrameMoving< double > &T) |
Apply a transformation (rotation and translation) represented by another ChFrameMoving T in local coordinate. More... | |
ChVector< double > | PointSpeedLocalToParent (const ChVector< double > &localpos) const |
Given the position of a point in local frame coords, and assuming it is sticky to frame, return the speed in parent coords. | |
ChVector< double > | PointSpeedLocalToParent (const ChVector< double > &localpos, const ChVector< double > &localspeed) const |
Given the position localpos of a point in the local reference frame, assuming that the point moves in the local reference frame with localspeed, return the speed in the parent reference frame. | |
ChVector< double > | PointAccelerationLocalToParent (const ChVector< double > &localpos) const |
Given the position of a point in local frame coords, and assuming it is sticky to frame, return the acceleration in parent coords. More... | |
ChVector< double > | PointAccelerationLocalToParent (const ChVector< double > &localpos, const ChVector< double > &localspeed, const ChVector< double > &localacc) const |
Given the position of a point in local frame coords, and assuming it has a frame-relative speed localspeed and frame-relative acceleration localacc, return the acceleration in parent coords. | |
ChVector< double > | PointSpeedParentToLocal (const ChVector< double > &parentpos, const ChVector< double > &parentspeed) const |
Given the position of a point in parent frame coords, and assuming it has an absolute speed parentspeed, return the speed in local coords. | |
ChVector< double > | PointAccelerationParentToLocal (const ChVector< double > &parentpos, const ChVector< double > &parentspeed, const ChVector< double > &parentacc) const |
Given the position of a point in parent frame coords, and assuming it has an absolute speed parentspeed and absolute acceleration parentacc, return the acceleration in local coords. | |
void | TransformLocalToParent (const ChFrameMoving< double > &local, ChFrameMoving< double > &parent) const |
This function transforms a frame from 'this' local coordinate system to parent frame coordinate system, and also transforms the speed and acceleration of the frame. More... | |
void | TransformParentToLocal (const ChFrameMoving< double > &parent, ChFrameMoving< double > &local) const |
This function transforms a frame from the parent coordinate system to 'this' local frame coordinate system. More... | |
bool | Equals (const ChFrameMoving< double > &other) const |
Returns true if coordsys is identical to other coordsys. | |
bool | Equals (const ChFrameMoving< double > &other, double tol) const |
Returns true if coordsys is equal to other coordsys, within a tolerance 'tol'. | |
virtual void | Invert () override |
The transformation (also for speeds, accelerations) is inverted in place. More... | |
ChFrameMoving< double > | GetInverse () const |
Public Member Functions inherited from chrono::ChFrame< double > | |
ChFrame (const ChVector< double > &mv=ChVector< double >(0, 0, 0), const ChQuaternion< double > &mq=ChQuaternion< double >(1, 0, 0, 0)) | |
Default constructor, or construct from pos and rot (as a quaternion) | |
ChFrame (const ChVector< double > &mv, const ChMatrix33< double > &ma) | |
Construct from pos and rotation (as a 3x3 matrix) | |
ChFrame (const ChCoordsys< double > &mc) | |
Construct from a coordsys. | |
ChFrame (const ChVector< double > &mv, const double alpha, const ChVector< double > &mu) | |
Construct from position mv and rotation of angle alpha around unit vector mu. | |
ChFrame (const ChFrame< double > &other) | |
Copy constructor, build from another frame. | |
ChFrame< double > & | operator= (const ChFrame< double > &other) |
Assignment operator: copy from another frame. | |
bool | operator== (const ChFrame< double > &other) const |
Returns true for identical frames. | |
bool | operator!= (const ChFrame< double > &other) const |
Returns true for different frames. | |
ChFrame< double > | operator>> (const ChFrame< double > &Fb) const |
The '>>' operator transforms a coordinate system, so transformations can be represented with this syntax: new_frame = old_frame >> tr_frame; For a sequence of transformations, i.e. More... | |
ChFrame< double > | operator* (const ChFrame< double > &Fb) const |
The '>>' operator transforms a vector, so transformations can be represented with this syntax: new_v = old_v >> tr_frame; For a sequence of transformations, i.e. More... | |
ChVector< double > | operator* (const ChVector< double > &V) const |
The '*' operator transforms a vector, so transformations can be represented with this syntax: new_v = tr_frame * old_v; For a sequence of transformations, i.e. More... | |
ChVector< double > | operator/ (const ChVector< double > &V) const |
The '/' is like the '*' operator (see), but uses the inverse transformation for A, in A/b. More... | |
ChFrame< double > & | operator>>= (const ChFrame< double > &T) |
Performs pre-multiplication of this frame by another frame, for example: A>>=T means A'=T*A ; or A'=A >> T. | |
ChFrame< double > & | operator>>= (const ChVector< double > &D) |
Performs pre-multiplication of this frame by a vector D, to 'move' by a displacement D: | |
ChFrame< double > & | operator>>= (const ChQuaternion< double > &R) |
Performs pre-multiplication of this frame by a quaternion R, to 'rotate' it by R: | |
ChFrame< double > & | operator>>= (const ChCoordsys< double > &F) |
Performs pre-multiplication of this frame by a ChCoordsys F, to transform it: | |
ChFrame< double > & | operator%= (const ChFrame< double > &T) |
Performs pre-multiplication of this frame by another frame, for example: A%=T means A'=T*A ; or A'=A >> T Note: DEPRECATED, use >>= instead. | |
ChFrame< double > & | operator*= (const ChFrame< double > &T) |
Performs post-multiplication of this frame by another frame, for example: A*=T means A'=A*T ; or A'=T >> A. | |
ChCoordsys< double > & | GetCoord () |
Return both current rotation and translation as a coordsystem object, with vector and quaternion. | |
const ChCoordsys< double > & | GetCoord () const |
ChVector< double > & | GetPos () |
Return the current translation as a 3d vector. | |
const ChVector< double > & | GetPos () const |
ChQuaternion< double > & | GetRot () |
Return the current rotation as a quaternion. | |
const ChQuaternion< double > & | GetRot () const |
ChMatrix33< double > & | GetA () |
Return the current rotation as a 3x3 matrix. | |
const ChMatrix33< double > & | GetA () const |
ChVector< double > | GetRotAxis () |
Get axis of finite rotation, in parent space. | |
double | GetRotAngle () |
Get angle of rotation about axis of finite rotation. | |
void | SetCoord (const ChCoordsys< double > &mcoord) |
Impose both translation and rotation as a single ChCoordsys. More... | |
void | SetCoord (const ChVector< double > &mv, const ChQuaternion< double > &mq) |
Impose both translation and rotation. More... | |
void | SetRot (const ChQuaternion< double > &mrot) |
Impose the rotation as a quaternion. More... | |
void | SetRot (const ChMatrix33< double > &mA) |
Impose the rotation as a 3x3 matrix. More... | |
void | SetPos (const ChVector< double > &mpos) |
Impose the translation. | |
void | ConcatenatePreTransformation (const ChFrame< double > &T) |
Apply a transformation (rotation and translation) represented by another ChFrame T. More... | |
void | ConcatenatePostTransformation (const ChFrame< double > &T) |
Apply a transformation (rotation and translation) represented by another ChFrame T in local coordinate. More... | |
void | Move (const ChVector< double > &V) |
An easy way to move the frame by the amount specified by vector V, (assuming V expressed in parent coordinates) | |
void | Move (const ChCoordsys< double > &VR) |
Apply both translation and rotation, assuming both expressed in parent coordinates, as a vector for translation and quaternion for rotation,. | |
ChVector< double > | TransformLocalToParent (const ChVector< double > &local) const |
This function transforms a point from the local frame coordinate system to the parent coordinate system. More... | |
void | TransformLocalToParent (const ChFrame< double > &local, ChFrame< double > &parent) const |
This function transforms a frame from 'this' local coordinate system to parent frame coordinate system. More... | |
ChVector< double > | TransformPointLocalToParent (const ChVector< double > &local) const |
ChVector< double > | TransformParentToLocal (const ChVector< double > &parent) const |
This function transforms a point from the parent coordinate system to local frame coordinate system. More... | |
void | TransformParentToLocal (const ChFrame< double > &parent, ChFrame< double > &local) const |
This function transforms a frame from the parent coordinate system to 'this' local frame coordinate system. More... | |
ChVector< double > | TransformPointParentToLocal (const ChVector< double > &parent) const |
ChVector< double > | TransformDirectionParentToLocal (const ChVector< double > &mdirection) const |
This function transforms a direction from 'this' local coordinate system to parent frame coordinate system. More... | |
ChVector< double > | TransformDirectionLocalToParent (const ChVector< double > &mdirection) const |
This function transforms a direction from the parent frame coordinate system to 'this' local coordinate system. More... | |
bool | Equals (const ChFrame< double > &other) const |
Returns true if coordsys is identical to other coordsys. | |
bool | Equals (const ChFrame< double > &other, double tol) const |
Returns true if coordsys is equal to other coordsys, within a tolerance 'tol'. | |
void | Normalize () |
Normalize the rotation, so that quaternion has unit length. | |
virtual void | SetIdentity () |
Sets to no translation and no rotation. | |
ChFrame< double > | GetInverse () const |
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. | |
Additional Inherited Members | |
Public Types inherited from chrono::ChVariableTupleCarrier_1vars< 6 > | |
typedef ChConstraintTuple_1vars< ChVariableTupleCarrier_1vars< N1 > > | type_constraint_tuple |
Public Attributes inherited from chrono::fea::ChNodeFEAbase | |
double | m_TotalMass |
Nodal mass obtained from element masss matrix. | |
Public Attributes inherited from chrono::ChFrameMoving< double > | |
ChCoordsys< double > | coord_dt |
Rotation and position speed, as vector+quaternion. | |
ChCoordsys< double > | coord_dtdt |
Rotation and position acceleration, as vector+quaternion. | |
Public Attributes inherited from chrono::ChFrame< double > | |
ChCoordsys< double > | coord |
Rotation and position, as vector+quaternion. | |
ChMatrix33< double > | Amatrix |
3x3 orthogonal rotation matrix | |
Static Public Attributes inherited from chrono::ChVariableTupleCarrier_1vars< 6 > | |
static const int | nvars1 |
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.
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, size is 6, it is {Force,Torque} both in absolute coords. 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
Implements chrono::ChLoadableUVW.
◆ 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
Implements chrono::ChLoadable.
◆ GetFixed()
|
inlineoverridevirtual |
Get the 'fixed' state of the node.
If true, its current field value is not changed by solver.
Implements chrono::fea::ChNodeFEAbase.
◆ SetFixed()
|
inlineoverridevirtual |
Set the 'fixed' state of the node.
If true, its current field value is not changed by solver.
Implements chrono::fea::ChNodeFEAbase.
◆ 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::ChNodeBase.
◆ 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::ChNodeBase.
◆ 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::ChNodeBase.
The documentation for this class was generated from the following files:
- /builds/uwsbel/chrono/src/chrono/fea/ChNodeFEAxyzrot.h
- /builds/uwsbel/chrono/src/chrono/fea/ChNodeFEAxyzrot.cpp