chrono::fea::ChElementTetraCorot_4_P Class Reference

Description

Tetrahedron FEM element with 4 nodes for scalar fields (for Poisson-like problems).

This is a classical element with linear displacement. EXPERIMENTAL

#include <ChElementTetraCorot_4.h>

Inheritance diagram for chrono::fea::ChElementTetraCorot_4_P:
Collaboration diagram for chrono::fea::ChElementTetraCorot_4_P:

Public Types

using ShapeVector = ChMatrixNM< double, 1, 4 >
 

Public Member Functions

virtual int GetNnodes () override
 Get the number of nodes used by this element.
 
virtual int GetNdofs () override
 Get the number of coordinates in the field used by the referenced nodes. More...
 
virtual int GetNodeNdofs (int n) override
 Get the number of coordinates from the specified node that are used by this element. More...
 
double GetVolume ()
 
virtual std::shared_ptr< ChNodeFEAbaseGetNodeN (int n) override
 Access the nth node.
 
virtual void SetNodes (std::shared_ptr< ChNodeFEAxyzP > nodeA, std::shared_ptr< ChNodeFEAxyzP > nodeB, std::shared_ptr< ChNodeFEAxyzP > nodeC, std::shared_ptr< ChNodeFEAxyzP > nodeD)
 
virtual void Update () override
 Update element at each time step.
 
virtual void ShapeFunctions (ShapeVector &N, double z0, double z1, double z2)
 Fills the N shape function matrix with the values of shape functions at zi parametric coordinates, where z0=1 at 1st vertex, z1=1 at second, z2 = 1 at third (volumetric shape functions). More...
 
virtual void GetStateBlock (ChVectorDynamic<> &mD) override
 Fills the D vector with the current field values at the nodes of the element, with proper ordering. More...
 
double ComputeVolume ()
 
virtual void ComputeStiffnessMatrix ()
 Computes the local STIFFNESS MATRIX of the element: K = Volume * [B]' * [D] * [B].
 
virtual void UpdateRotation () override
 Given the actual position of the nodes, recompute the cumulative rotation matrix A.
 
virtual void ComputeKRMmatricesGlobal (ChMatrixRef H, double Kfactor, double Rfactor=0, double Mfactor=0) override
 Sets H as the global stiffness matrix K, scaled by Kfactor. More...
 
virtual void ComputeInternalForces (ChVectorDynamic<> &Fi) override
 Computes the internal 'pseudo-forces' and set values in the Fi vector.
 
void SetMaterial (std::shared_ptr< ChContinuumPoisson3D > my_material)
 Set the material of the element.
 
std::shared_ptr< ChContinuumPoisson3DGetMaterial ()
 
const ChMatrixDynamicGetMatrB () const
 Get the partial derivatives matrix MatrB and the StiffnessMatrix.
 
const ChMatrixDynamicGetStiffnessMatrix () const
 
ChVectorN< double, 3 > GetPgradient ()
 Returns the gradient of P (note that the tetrahedron 4 nodes is a linear element, thus the gradient is constant in the entire volume). 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: here the {t} temperature.
 
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 N'*F , where N is some type of shape function evaluated at U,V,W coordinates of the volume, each ranging in 0..+1 F is a load, N'*F is the resulting generalized load Returns also det[J] with J=[dx/du,..], that might be useful in gauss quadrature. More...
 
virtual double GetDensity () override
 Return 0 if not supportable by ChLoaderVolumeGravity.
 
virtual bool IsTetrahedronIntegrationNeeded () override
 If true, use quadrature over u,v,w in [0..1] range as tetrahedron volumetric coords, with z=1-u-v-w otherwise use quadrature over u,v,w in [-1..+1] as box isoparametric coords.
 
- Public Member Functions inherited from chrono::fea::ChElementGeneric
ChKblockGenericKstiffness ()
 Access the proxy to stiffness, for sparse solver.
 
virtual void EleIntLoadResidual_F (ChVectorDynamic<> &R, const double c) override
 Add the internal forces (pasted at global nodes offsets) into a global vector R, multiplied by a scaling factor c, as R += forces * c This default implementation is SLIGHTLY INEFFICIENT.
 
virtual void EleIntLoadResidual_Mv (ChVectorDynamic<> &R, const ChVectorDynamic<> &w, const double c) override
 Add the product of element mass M by a vector w (pasted at global nodes offsets) into a global vector R, multiplied by a scaling factor c, as R += M * w * c This default implementation is VERY INEFFICIENT.
 
virtual void EleIntLoadResidual_F_gravity (ChVectorDynamic<> &R, const ChVector<> &G_acc, const double c) override
 Add the contribution of gravity loads, multiplied by a scaling factor c, as: R += M * g * c This default implementation is VERY INEFFICIENT. More...
 
virtual void ComputeGravityForces (ChVectorDynamic<> &Fg, const ChVector<> &G_acc) override
 Compute the gravitational forces. More...
 
virtual void ComputeMmatrixGlobal (ChMatrixRef M) override
 Calculate the mass matrix, expressed in global reference. More...
 
virtual void InjectKRMmatrices (ChSystemDescriptor &descriptor) override
 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) override
 Add the current stiffness K and damping R and mass M matrices in encapsulated ChKblock item(s), if any. More...
 
virtual void VariablesFbLoadInternalForces (double factor=1.) override
 Add the internal forces, expressed as nodal forces, into the encapsulated ChVariables.
 
virtual void VariablesFbIncrementMq () override
 Add M*q (internal masses multiplied current 'qb').
 
- Public Member Functions inherited from chrono::fea::ChElementBase
virtual int GetNdofs_active ()
 Get the actual number of active degrees of freedom. More...
 
virtual int GetNodeNdofs_active (int n)
 Get the actual number of active coordinates from the specified node that are used by this element. More...
 
virtual void ComputeNodalMass ()
 Compute element's nodal masses.
 
virtual void EleDoIntegration ()
 This is optionally implemented if there is some internal state that requires integration.
 
- Public Member Functions inherited from chrono::fea::ChElementCorotational
ChMatrix33Rotation ()
 Access the cumulative rotation matrix of the element. More...
 
- Public Member Functions inherited from chrono::ChLoadableUVW
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

- Protected Attributes inherited from chrono::fea::ChElementGeneric
ChKblockGeneric Kmatr
 
- Protected Attributes inherited from chrono::fea::ChElementCorotational
ChMatrix33 A
 

Member Function Documentation

◆ ComputeKRMmatricesGlobal()

void chrono::fea::ChElementTetraCorot_4_P::ComputeKRMmatricesGlobal ( ChMatrixRef  H,
double  Kfactor,
double  Rfactor = 0,
double  Mfactor = 0 
)
overridevirtual

Sets H as the global stiffness matrix K, scaled by Kfactor.

Optionally, also superimposes global damping matrix R, scaled by Rfactor, and global mass matrix M multiplied by Mfactor.

Implements chrono::fea::ChElementBase.

◆ ComputeNF()

void chrono::fea::ChElementTetraCorot_4_P::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 N'*F , where N is some type of shape function evaluated at U,V,W coordinates of the volume, each ranging in 0..+1 F is a load, N'*F is the resulting generalized load Returns also det[J] with J=[dx/du,..], that might be useful in gauss quadrature.

Parameters
Uparametric coordinate in volume
Vparametric coordinate in volume
Wparametric coordinate in volume
QiReturn result of N'*F here, maybe with offset block_offset
detJReturn det[J] here
FInput F vector, size is = n.field coords.
state_xif != 0, update state (pos. part) to this, then evaluate Q
state_wif != 0, update state (speed part) to this, then evaluate Q

Implements chrono::ChLoadableUVW.

◆ GetNdofs()

virtual int chrono::fea::ChElementTetraCorot_4_P::GetNdofs ( )
inlineoverridevirtual

Get the number of coordinates in the field used by the referenced nodes.

This is for example the size (number of rows/columns) of the local stiffness matrix.

Implements chrono::fea::ChElementBase.

◆ GetNodeNdofs()

virtual int chrono::fea::ChElementTetraCorot_4_P::GetNodeNdofs ( int  n)
inlineoverridevirtual

Get the number of coordinates from the specified node that are used by this element.

Note that this may be different from the value returned by GetNodeN(n)->GetNdofW().

Implements chrono::fea::ChElementBase.

◆ GetPgradient()

ChVectorN< double, 3 > chrono::fea::ChElementTetraCorot_4_P::GetPgradient ( )

Returns the gradient of P (note that the tetrahedron 4 nodes is a linear element, thus the gradient is constant in the entire volume).

It is in the original undeformed unrotated reference.

◆ GetStateBlock()

void chrono::fea::ChElementTetraCorot_4_P::GetStateBlock ( ChVectorDynamic<> &  mD)
overridevirtual

Fills the D vector with the current field values at the nodes of the element, with proper ordering.

If the D vector has not the size of this->GetNdofs(), it will be resized. For corotational elements, field is assumed in local reference!

Implements chrono::fea::ChElementBase.

◆ ShapeFunctions()

void chrono::fea::ChElementTetraCorot_4_P::ShapeFunctions ( ShapeVector &  N,
double  z0,
double  z1,
double  z2 
)
virtual

Fills the N shape function matrix with the values of shape functions at zi parametric coordinates, where z0=1 at 1st vertex, z1=1 at second, z2 = 1 at third (volumetric shape functions).

The 4th is computed form the first 3. All ranging in [0...1]. NOTE! actually N should be a 3row, 12 column sparse matrix, as N = [n1*eye(3) n2*eye(3) n3*eye(3) n4*eye(3)]; , but to avoid wasting zero and repeated elements, here it stores only the n1 n2 n3 n4 values in a 1 row, 4 columns matrix!


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