KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
Public Types | Public Member Functions | Static Protected Attributes | Friends | List of all members
Kratos::CrBeamElement3D2N Class Reference

This is a 3D-2node beam element with 3 translational dofs and 3 rotational dof per node. More...

#include <cr_beam_element_3D2N.hpp>

Inheritance diagram for Kratos::CrBeamElement3D2N:
Collaboration diagram for Kratos::CrBeamElement3D2N:

Public Types

typedef Element BaseType
 
typedef BaseType::GeometryType GeometryType
 
typedef BaseType::NodesArrayType NodesArrayType
 
typedef BaseType::PropertiesType PropertiesType
 
typedef BaseType::IndexType IndexType
 
typedef BaseType::SizeType SizeType
 
typedef BaseType::MatrixType MatrixType
 
typedef BaseType::VectorType VectorType
 
typedef BaseType::EquationIdVectorType EquationIdVectorType
 
typedef BaseType::DofsVectorType DofsVectorType
 
- Public Types inherited from Kratos::Element
typedef Element ElementType
 definition of element type More...
 
typedef GeometricalObject BaseType
 base type: an GeometricalObject that automatically has a unique number More...
 
typedef Node NodeType
 definition of node type (default is: Node) More...
 
typedef Properties PropertiesType
 
typedef Geometry< NodeTypeGeometryType
 definition of the geometry type with given NodeType More...
 
typedef Geometry< NodeType >::PointsArrayType NodesArrayType
 definition of nodes container type, redefined from GeometryType More...
 
typedef Vector VectorType
 
typedef Matrix MatrixType
 
typedef std::size_t IndexType
 
typedef std::size_t SizeType
 
typedef Dof< doubleDofType
 
typedef std::vector< std::size_t > EquationIdVectorType
 
typedef std::vector< DofType::PointerDofsVectorType
 
typedef PointerVectorSet< DofTypeDofsArrayType
 
typedef GeometryData::IntegrationMethod IntegrationMethod
 Type definition for integration methods. More...
 
typedef GeometryData GeometryDataType
 
- Public Types inherited from Kratos::GeometricalObject
typedef Node NodeType
 Definition of the node type. More...
 
typedef Geometry< NodeTypeGeometryType
 The geometry type definition. More...
 
typedef std::size_t IndexType
 Defines the index type. More...
 
typedef std::size_t result_type
 Defines the result type. More...
 
- Public Types inherited from Kratos::IndexedObject
typedef std::size_t IndexType
 The definition of the index type. More...
 
typedef std::size_t result_type
 The definition of the result_type. More...
 
- Public Types inherited from Kratos::Flags
enum  FlagsList {
  Flag0 = BlockType(1) , Flag1 = BlockType(1) << 1 , Flag2 = BlockType(1) << 2 , Flag3 = BlockType(1) << 3 ,
  Flag4 = BlockType(1) << 4 , Flag5 = BlockType(1) << 5 , Flag6 = BlockType(1) << 6 , Flag7 = BlockType(1) << 7 ,
  Flag8 = BlockType(1) << 8 , Flag9 = BlockType(1) << 9 , Flag10 = BlockType(1) << 10 , Flag11 = BlockType(1) << 11 ,
  Flag12 = BlockType(1) << 12 , Flag13 = BlockType(1) << 13 , Flag14 = BlockType(1) << 14 , Flag15 = BlockType(1) << 15 ,
  Flag16 = BlockType(1) << 16 , Flag17 = BlockType(1) << 17 , Flag18 = BlockType(1) << 18 , Flag19 = BlockType(1) << 19 ,
  Flag20 = BlockType(1) << 20 , Flag21 = BlockType(1) << 21 , Flag22 = BlockType(1) << 22 , Flag23 = BlockType(1) << 23 ,
  Flag24 = BlockType(1) << 24 , Flag25 = BlockType(1) << 25 , Flag26 = BlockType(1) << 26 , Flag27 = BlockType(1) << 27 ,
  Flag28 = BlockType(1) << 28 , Flag29 = BlockType(1) << 29 , Flag30 = BlockType(1) << 30
}
 
typedef int64_t BlockType
 
typedef int64_t FlagType
 
typedef std::size_t IndexType
 

Public Member Functions

 KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION (CrBeamElement3D2N)
 
 CrBeamElement3D2N ()
 
 CrBeamElement3D2N (IndexType NewId, GeometryType::Pointer pGeometry)
 
 CrBeamElement3D2N (IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
 
 ~CrBeamElement3D2N () override
 
Element::Pointer Create (IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
 Creates a new element. More...
 
Element::Pointer Create (IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
 Creates a new element. More...
 
void EquationIdVector (EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
 
void GetDofList (DofsVectorType &rElementalDofList, const ProcessInfo &rCurrentProcessInfo) const override
 
BoundedMatrix< double, msElementSize, msElementSizeCreateElementStiffnessMatrix_Material () const
 This function calculates the elastic part of the total stiffness matrix. More...
 
BoundedMatrix< double, msElementSize, msElementSizeCreateElementStiffnessMatrix_Geometry () const
 This function calculates the geometric part of the total stiffness matrix. More...
 
virtual BoundedMatrix< double, msLocalSize, msLocalSizeCalculateDeformationStiffness () const
 This function calculates the element stiffness w.r.t. deformation modes. More...
 
BoundedMatrix< double, msElementSize, msLocalSizeCalculateTransformationS () const
 This function calculates a transformation matrix from deformation modes to real deformations. More...
 
BoundedVector< double, msLocalSizeGetCurrentNodalPosition () const
 This function calculates the current nodal position. More...
 
BoundedVector< double, msLocalSizeCalculateElementForces () const
 This function calculates the internal element forces. More...
 
BoundedMatrix< double, msElementSize, msElementSizeCalculateInitialLocalCS () const
 This function calculates the initial transformation matrix to globalize/localize vectors and/or matrices. More...
 
BoundedMatrix< double, msDimension, msDimensionUpdateRotationMatrixLocal (Vector &Bisectrix, Vector &VectorDifference) const
 This function updates constantly the transformation matrix. More...
 
void SaveQuaternionParameters ()
 
void CalculateLocalSystem (MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
 
void ConstCalculateLocalSystem (MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) const
 
void CalculateRightHandSide (VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
 
void CalculateLeftHandSide (MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
 
virtual void ConstCalculateRightHandSide (VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) const
 
void ConstCalculateLeftHandSide (MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) const
 
void CalculateMassMatrix (MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) override
 
void CalculateLumpedMassMatrix (MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) const
 This function calculates the lumped mass matrix. More...
 
void CalculateConsistentMassMatrix (MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) const
 This function calculates the consistent mass matrix. More...
 
void BuildSingleMassMatrix (MatrixType &rMassMatrix, const double Phi, const double CT, const double CR, const double L, const double dir) const
 This function calculates parts of the total consistent mass matrix to simplify the code. More...
 
void CalculateDampingMatrix (MatrixType &rDampingMatrix, const ProcessInfo &rCurrentProcessInfo) override
 
void AddExplicitContribution (const VectorType &rRHSVector, const Variable< VectorType > &rRHSVariable, const Variable< array_1d< double, 3 > > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo) override
 This function is designed to make the element to assemble an rRHS vector identified by a variable rRHSVariable by assembling it to the nodes on the variable rDestinationVariable. (This is the vector version) More...
 
void GetValuesVector (Vector &rValues, int Step=0) const override
 
void GetSecondDerivativesVector (Vector &rValues, int Step=0) const override
 
void GetFirstDerivativesVector (Vector &rValues, int Step=0) const override
 
void AssembleSmallInBigMatrix (const Matrix &rSmallMatrix, BoundedMatrix< double, msElementSize, msElementSize > &rBigMatrix) const
 This function is used to assemble single transformation matrix in the big global rotation matrix. More...
 
int Check (const ProcessInfo &rCurrentProcessInfo) const override
 
double CalculatePsi (const double I, const double A_eff) const
 This function calculates reduction values in case of shear-deformable structures. More...
 
double CalculateShearModulus () const
 This function calculates shear modulus from user input values. More...
 
BoundedVector< double, msElementSizeCalculateBodyForces () const
 This function calculates self-weight forces. More...
 
void Calculate (const Variable< Matrix > &rVariable, Matrix &rOutput, const ProcessInfo &rCurrentProcessInfo) override
 
void CalculateOnIntegrationPoints (const Variable< array_1d< double, 3 > > &rVariable, std::vector< array_1d< double, 3 > > &rOutput, const ProcessInfo &rCurrentProcessInfo) override
 
IntegrationMethod GetIntegrationMethod () const override
 
void CalculateAndAddWorkEquivalentNodalForcesLineLoad (const BoundedVector< double, msDimension > ForceInput, BoundedVector< double, msElementSize > &rRightHandSideVector, const double GeometryLength) const
 This function calculates nodal moments due to self-weight. More...
 
Vector CalculateSymmetricDeformationMode () const
 This function calculates the symmetric deformation modes. More...
 
Vector CalculateAntiSymmetricDeformationMode () const
 This function calculates the antisymmetric deformation modes. More...
 
Vector CalculateLocalNodalForces () const
 This function calculates the local nodal forces. More...
 
void FinalizeNonLinearIteration (const ProcessInfo &rCurrentProcessInfo) override
 
Vector CalculateGlobalNodalForces () const
 
Vector GetIncrementDeformation () const
 
BoundedMatrix< double, msElementSize, msElementSizeGetTransformationMatrixGlobal () const
 
void InitializeNonLinearIteration (const ProcessInfo &rCurrentProcessInfo) override
 
void UpdateQuaternionParameters (double &rScalNodeA, double &rScalNodeB, Vector &rVecNodeA, Vector &rVecNodeB) const
 
const Parameters GetSpecifications () const override
 This method provides the specifications/requirements of the element. More...
 
- Public Member Functions inherited from Kratos::Element
 Element (IndexType NewId=0)
 
 Element (IndexType NewId, const NodesArrayType &ThisNodes)
 
 Element (IndexType NewId, GeometryType::Pointer pGeometry)
 
 Element (IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
 
 Element (Element const &rOther)
 Copy constructor. More...
 
 ~Element () override
 Destructor. More...
 
Elementoperator= (Element const &rOther)
 Assignment operator. More...
 
virtual Pointer Clone (IndexType NewId, NodesArrayType const &ThisNodes) const
 It creates a new element pointer and clones the previous element data. More...
 
virtual void Initialize (const ProcessInfo &rCurrentProcessInfo)
 
virtual void ResetConstitutiveLaw ()
 
virtual void InitializeSolutionStep (const ProcessInfo &rCurrentProcessInfo)
 
virtual void FinalizeSolutionStep (const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateFirstDerivativesContributions (MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateFirstDerivativesLHS (MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateFirstDerivativesRHS (VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateSecondDerivativesContributions (MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateSecondDerivativesLHS (MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateSecondDerivativesRHS (VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateLumpedMassVector (VectorType &rLumpedMassVector, const ProcessInfo &rCurrentProcessInfo) const
 
virtual void AddExplicitContribution (const ProcessInfo &rCurrentProcessInfo)
 
virtual void AddExplicitContribution (const VectorType &rRHSVector, const Variable< VectorType > &rRHSVariable, const Variable< double > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo)
 This function is designed to make the element to assemble an rRHS vector identified by a variable rRHSVariable by assembling it to the nodes on the variable rDestinationVariable. (This is the double version) More...
 
virtual void AddExplicitContribution (const MatrixType &rLHSMatrix, const Variable< MatrixType > &rLHSVariable, const Variable< Matrix > &rDestinationVariable, const ProcessInfo &rCurrentProcessInfo)
 This function is designed to make the element to assemble an rRHS vector identified by a variable rRHSVariable by assembling it to the nodes on the variable rDestinationVariable. (This is the matrix version) More...
 
virtual void Calculate (const Variable< double > &rVariable, double &Output, const ProcessInfo &rCurrentProcessInfo)
 
virtual void Calculate (const Variable< array_1d< double, 3 > > &rVariable, array_1d< double, 3 > &Output, const ProcessInfo &rCurrentProcessInfo)
 
virtual void Calculate (const Variable< Vector > &rVariable, Vector &Output, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateOnIntegrationPoints (const Variable< bool > &rVariable, std::vector< bool > &rOutput, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateOnIntegrationPoints (const Variable< int > &rVariable, std::vector< int > &rOutput, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateOnIntegrationPoints (const Variable< double > &rVariable, std::vector< double > &rOutput, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateOnIntegrationPoints (const Variable< array_1d< double, 3 >> &rVariable, std::vector< array_1d< double, 3 >> &rOutput, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateOnIntegrationPoints (const Variable< array_1d< double, 4 >> &rVariable, std::vector< array_1d< double, 4 >> &rOutput, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateOnIntegrationPoints (const Variable< array_1d< double, 6 >> &rVariable, std::vector< array_1d< double, 6 >> &rOutput, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateOnIntegrationPoints (const Variable< array_1d< double, 9 >> &rVariable, std::vector< array_1d< double, 9 >> &rOutput, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateOnIntegrationPoints (const Variable< Vector > &rVariable, std::vector< Vector > &rOutput, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateOnIntegrationPoints (const Variable< Matrix > &rVariable, std::vector< Matrix > &rOutput, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateOnIntegrationPoints (const Variable< ConstitutiveLaw::Pointer > &rVariable, std::vector< ConstitutiveLaw::Pointer > &rOutput, const ProcessInfo &rCurrentProcessInfo)
 
virtual void SetValuesOnIntegrationPoints (const Variable< bool > &rVariable, const std::vector< bool > &rValues, const ProcessInfo &rCurrentProcessInfo)
 
virtual void SetValuesOnIntegrationPoints (const Variable< int > &rVariable, const std::vector< int > &rValues, const ProcessInfo &rCurrentProcessInfo)
 
virtual void SetValuesOnIntegrationPoints (const Variable< double > &rVariable, const std::vector< double > &rValues, const ProcessInfo &rCurrentProcessInfo)
 
virtual void SetValuesOnIntegrationPoints (const Variable< array_1d< double, 3 >> &rVariable, const std::vector< array_1d< double, 3 >> &rValues, const ProcessInfo &rCurrentProcessInfo)
 
virtual void SetValuesOnIntegrationPoints (const Variable< array_1d< double, 4 >> &rVariable, const std::vector< array_1d< double, 4 >> &rValues, const ProcessInfo &rCurrentProcessInfo)
 
virtual void SetValuesOnIntegrationPoints (const Variable< array_1d< double, 6 >> &rVariable, const std::vector< array_1d< double, 6 >> &rValues, const ProcessInfo &rCurrentProcessInfo)
 
virtual void SetValuesOnIntegrationPoints (const Variable< array_1d< double, 9 >> &rVariable, const std::vector< array_1d< double, 9 >> &rValues, const ProcessInfo &rCurrentProcessInfo)
 
virtual void SetValuesOnIntegrationPoints (const Variable< Vector > &rVariable, const std::vector< Vector > &rValues, const ProcessInfo &rCurrentProcessInfo)
 
virtual void SetValuesOnIntegrationPoints (const Variable< Matrix > &rVariable, const std::vector< Matrix > &rValues, const ProcessInfo &rCurrentProcessInfo)
 
virtual void SetValuesOnIntegrationPoints (const Variable< ConstitutiveLaw::Pointer > &rVariable, const std::vector< ConstitutiveLaw::Pointer > &rValues, const ProcessInfo &rCurrentProcessInfo)
 
virtual void MassMatrix (MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo)
 
virtual void AddMassMatrix (MatrixType &rLeftHandSideMatrix, double coeff, const ProcessInfo &rCurrentProcessInfo)
 
virtual void DampMatrix (MatrixType &rDampMatrix, const ProcessInfo &rCurrentProcessInfo)
 
virtual void AddInertiaForces (VectorType &rRightHandSideVector, double coeff, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateLocalVelocityContribution (MatrixType &rDampingMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateSensitivityMatrix (const Variable< double > &rDesignVariable, Matrix &rOutput, const ProcessInfo &rCurrentProcessInfo)
 
virtual void CalculateSensitivityMatrix (const Variable< array_1d< double, 3 > > &rDesignVariable, Matrix &rOutput, const ProcessInfo &rCurrentProcessInfo)
 
PropertiesType::Pointer pGetProperties ()
 returns the pointer to the property of the element. Does not throw an error, to allow copying of elements which don't have any property assigned. More...
 
const PropertiesType::Pointer pGetProperties () const
 
PropertiesTypeGetProperties ()
 
PropertiesType const & GetProperties () const
 
void SetProperties (PropertiesType::Pointer pProperties)
 
bool HasProperties () const
 Check that the Element has a correctly initialized pointer to a Properties instance. More...
 
std::string Info () const override
 Turn back information as a string. More...
 
void PrintInfo (std::ostream &rOStream) const override
 Print information about this object. More...
 
void PrintData (std::ostream &rOStream) const override
 Print object's data. More...
 
 KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION (Element)
 
- Public Member Functions inherited from Kratos::GeometricalObject
 GeometricalObject (IndexType NewId=0)
 Default constructor. More...
 
 GeometricalObject (IndexType NewId, GeometryType::Pointer pGeometry)
 Default constructor. More...
 
 ~GeometricalObject () override
 Destructor. More...
 
 GeometricalObject (GeometricalObject const &rOther)
 Copy constructor. More...
 
GeometricalObjectoperator= (GeometricalObject const &rOther)
 Assignment operator. More...
 
virtual void SetGeometry (GeometryType::Pointer pGeometry)
 Sets the pointer to the geometry. More...
 
GeometryType::Pointer pGetGeometry ()
 Returns the pointer to the geometry. More...
 
const GeometryType::Pointer pGetGeometry () const
 Returns the pointer to the geometry (const version) More...
 
GeometryTypeGetGeometry ()
 Returns the reference of the geometry. More...
 
GeometryType const & GetGeometry () const
 Returns the reference of the geometry (const version) More...
 
FlagsGetFlags ()
 Returns the flags of the object. More...
 
Flags const & GetFlags () const
 Returns the flags of the object (const version) More...
 
void SetFlags (Flags const &rThisFlags)
 Sets the flags of the object. More...
 
DataValueContainerData ()
 
DataValueContainerGetData ()
 
DataValueContainer const & GetData () const
 
void SetData (DataValueContainer const &rThisData)
 
template<class TDataType >
bool Has (const Variable< TDataType > &rThisVariable) const
 
template<class TVariableType >
void SetValue (const TVariableType &rThisVariable, typename TVariableType::Type const &rValue)
 
template<class TVariableType >
TVariableType::Type & GetValue (const TVariableType &rThisVariable)
 
template<class TVariableType >
TVariableType::Type const & GetValue (const TVariableType &rThisVariable) const
 
std::string Info () const override
 Turn back information as a string. More...
 
void PrintInfo (std::ostream &rOStream) const override
 Print information about this object. More...
 
void PrintData (std::ostream &rOStream) const override
 Print object's data. More...
 
unsigned int use_count () const noexcept
 
 KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION (GeometricalObject)
 Pointer definition of GeometricalObject. More...
 
bool IsActive () const
 Checks if the GeometricalObject is active. More...
 
- Public Member Functions inherited from Kratos::IndexedObject
 IndexedObject (IndexType NewId=0)
 Default constructor. More...
 
virtual ~IndexedObject ()
 Destructor. More...
 
 IndexedObject (IndexedObject const &rOther)
 Copy constructor. More...
 
IndexedObjectoperator= (IndexedObject const &rOther)
 Assignment operator. More...
 
template<class TObjectType >
IndexType operator() (TObjectType const &rThisObject) const
 
IndexType Id () const
 
IndexType GetId () const
 
virtual void SetId (IndexType NewId)
 
IndexTypeDepricatedIdAccess ()
 TODO: remove this function when removing data_file_io object. More...
 
 KRATOS_CLASS_POINTER_DEFINITION (IndexedObject)
 Pointer definition of IndexedObject. More...
 
- Public Member Functions inherited from Kratos::Flags
Flagsoperator= (Flags const &rOther)
 Assignment operator. More...
 
 operator bool () const
 
Flags operator~ () const
 
bool operator! () const
 
void AssignFlags (Flags const &rOther)
 
void Set (const Flags ThisFlag)
 
void Set (const Flags ThisFlag, bool Value)
 
void Reset (const Flags ThisFlag)
 
void Flip (const Flags ThisFlag)
 
void SetPosition (IndexType Position, bool Value=true)
 
bool GetPosition (IndexType Position) const
 
void FlipPosition (IndexType Position)
 
void ClearPosition (IndexType Position)
 
void Clear ()
 
Flags AsFalse () const
 
bool Is (Flags const &rOther) const
 
bool IsDefined (Flags const &rOther) const
 
bool IsNot (Flags const &rOther) const
 
bool IsNotDefined (Flags const &rOther) const
 
 KRATOS_CLASS_POINTER_DEFINITION (Flags)
 Pointer definition of Flags. More...
 
const Flagsoperator|= (const Flags &Other)
 
const Flagsoperator&= (const Flags &Other)
 
 Flags ()
 Default constructor. More...
 
 Flags (Flags const &rOther)
 Copy constructor. More...
 
virtual ~Flags ()
 Destructor. More...
 

Static Protected Attributes

static constexpr int msNumberOfNodes = 2
 
static constexpr int msDimension = 3
 
static constexpr unsigned int msLocalSize = msNumberOfNodes * msDimension
 
static constexpr unsigned int msElementSize = msLocalSize * 2
 

Friends

class Serializer
 

Additional Inherited Members

- Static Public Member Functions inherited from Kratos::GeometricalObject
static bool HasSameType (const GeometricalObject &rLHS, const GeometricalObject &rRHS)
 Checks if two GeometricalObject have the same type. More...
 
static bool HasSameType (const GeometricalObject *rLHS, const GeometricalObject *rRHS)
 Checks if two GeometricalObject have the same type (pointer version) More...
 
static bool HasSameGeometryType (const GeometricalObject &rLHS, const GeometricalObject &rRHS)
 Checks if two GeometricalObject have the same geometry type. More...
 
static bool HasSameGeometryType (const GeometricalObject *rLHS, const GeometricalObject *rRHS)
 Checks if two GeometricalObject have the same geometry type (pointer version) More...
 
static bool IsSame (const GeometricalObject &rLHS, const GeometricalObject &rRHS)
 Checks if two GeometricalObject are the same. More...
 
static bool IsSame (const GeometricalObject *rLHS, const GeometricalObject *rRHS)
 Checks if two GeometricalObject are the same (pointer version) More...
 
- Static Public Member Functions inherited from Kratos::Flags
static const Flags AllDefined ()
 
static const Flags AllTrue ()
 
static Flags Create (IndexType ThisPosition, bool Value=true)
 

Detailed Description

This is a 3D-2node beam element with 3 translational dofs and 3 rotational dof per node.

Author
Klaus B Sautter

Member Typedef Documentation

◆ BaseType

◆ DofsVectorType

typedef BaseType::DofsVectorType Kratos::CrBeamElement3D2N::DofsVectorType

◆ EquationIdVectorType

typedef BaseType::EquationIdVectorType Kratos::CrBeamElement3D2N::EquationIdVectorType

◆ GeometryType

◆ IndexType

◆ MatrixType

typedef BaseType::MatrixType Kratos::CrBeamElement3D2N::MatrixType

◆ NodesArrayType

typedef BaseType::NodesArrayType Kratos::CrBeamElement3D2N::NodesArrayType

◆ PropertiesType

typedef BaseType::PropertiesType Kratos::CrBeamElement3D2N::PropertiesType

◆ SizeType

typedef BaseType::SizeType Kratos::CrBeamElement3D2N::SizeType

◆ VectorType

typedef BaseType::VectorType Kratos::CrBeamElement3D2N::VectorType

Constructor & Destructor Documentation

◆ CrBeamElement3D2N() [1/3]

Kratos::CrBeamElement3D2N::CrBeamElement3D2N ( )
inline

◆ CrBeamElement3D2N() [2/3]

Kratos::CrBeamElement3D2N::CrBeamElement3D2N ( IndexType  NewId,
GeometryType::Pointer  pGeometry 
)

◆ CrBeamElement3D2N() [3/3]

Kratos::CrBeamElement3D2N::CrBeamElement3D2N ( IndexType  NewId,
GeometryType::Pointer  pGeometry,
PropertiesType::Pointer  pProperties 
)

◆ ~CrBeamElement3D2N()

Kratos::CrBeamElement3D2N::~CrBeamElement3D2N ( )
override

Member Function Documentation

◆ AddExplicitContribution()

void Kratos::CrBeamElement3D2N::AddExplicitContribution ( const VectorType rRHSVector,
const Variable< VectorType > &  rRHSVariable,
const Variable< array_1d< double, 3 > > &  rDestinationVariable,
const ProcessInfo rCurrentProcessInfo 
)
overridevirtual

This function is designed to make the element to assemble an rRHS vector identified by a variable rRHSVariable by assembling it to the nodes on the variable rDestinationVariable. (This is the vector version)

The "AddExplicit" FUNCTIONS THE ONLY FUNCTIONS IN WHICH AN ELEMENT IS ALLOWED TO WRITE ON ITS NODES. The caller is expected to ensure thread safety hence SET-/UNSET-LOCK MUST BE PERFORMED IN THE STRATEGY BEFORE CALLING THIS FUNCTION

Parameters
rRHSVectorinput variable containing the RHS vector to be assembled
rRHSVariablevariable describing the type of the RHS vector to be assembled
rDestinationVariablevariable in the database to which the rRHSvector will be assembled
rCurrentProcessInfothe current process info instance

Reimplemented from Kratos::Element.

◆ AssembleSmallInBigMatrix()

void Kratos::CrBeamElement3D2N::AssembleSmallInBigMatrix ( const Matrix rSmallMatrix,
BoundedMatrix< double, msElementSize, msElementSize > &  rBigMatrix 
) const

This function is used to assemble single transformation matrix in the big global rotation matrix.

Parameters
rSmallMatrixThe local transformation matrix
rBigMatrixThe total global rotation matrix

◆ BuildSingleMassMatrix()

void Kratos::CrBeamElement3D2N::BuildSingleMassMatrix ( MatrixType rMassMatrix,
const double  Phi,
const double  CT,
const double  CR,
const double  L,
const double  dir 
) const

This function calculates parts of the total consistent mass matrix to simplify the code.

Parameters
rMassMatrixThe current mass matrix
PhiThe reduction value in case of shear-deformable structures
CTA scaling factor
CRA scaling factor
LThe element length
dirThe direction of the current cs

◆ Calculate()

void Kratos::CrBeamElement3D2N::Calculate ( const Variable< Matrix > &  rVariable,
Matrix rOutput,
const ProcessInfo rCurrentProcessInfo 
)
overridevirtual

Reimplemented from Kratos::Element.

Reimplemented in Kratos::CrBeamElementLinear3D2N.

◆ CalculateAndAddWorkEquivalentNodalForcesLineLoad()

void Kratos::CrBeamElement3D2N::CalculateAndAddWorkEquivalentNodalForcesLineLoad ( const BoundedVector< double, msDimension ForceInput,
BoundedVector< double, msElementSize > &  rRightHandSideVector,
const double  GeometryLength 
) const

This function calculates nodal moments due to self-weight.

Parameters
ForceInputThe self-weight line load vector
rRightHandSideVectorThe right hand side of the problem
GeometryLengthThe element length

◆ CalculateAntiSymmetricDeformationMode()

Vector Kratos::CrBeamElement3D2N::CalculateAntiSymmetricDeformationMode ( ) const

This function calculates the antisymmetric deformation modes.

◆ CalculateBodyForces()

BoundedVector< double, CrBeamElement3D2N::msElementSize > Kratos::CrBeamElement3D2N::CalculateBodyForces ( ) const

This function calculates self-weight forces.

◆ CalculateConsistentMassMatrix()

void Kratos::CrBeamElement3D2N::CalculateConsistentMassMatrix ( MatrixType rMassMatrix,
const ProcessInfo rCurrentProcessInfo 
) const

This function calculates the consistent mass matrix.

Parameters
rMassMatrixThe current mass matrix
rCurrentProcessInfoThe current Process information

◆ CalculateDampingMatrix()

void Kratos::CrBeamElement3D2N::CalculateDampingMatrix ( MatrixType rDampingMatrix,
const ProcessInfo rCurrentProcessInfo 
)
overridevirtual

this is called during the assembling process in order to calculate the elemental damping matrix

Parameters
rDampingMatrixthe elemental damping matrix
rCurrentProcessInfothe current process info instance

Reimplemented from Kratos::Element.

◆ CalculateDeformationStiffness()

BoundedMatrix< double, CrBeamElement3D2N::msLocalSize, CrBeamElement3D2N::msLocalSize > Kratos::CrBeamElement3D2N::CalculateDeformationStiffness ( ) const
virtual

This function calculates the element stiffness w.r.t. deformation modes.

Reimplemented in Kratos::CrBeamElementLinear3D2N.

◆ CalculateElementForces()

BoundedVector< double, CrBeamElement3D2N::msLocalSize > Kratos::CrBeamElement3D2N::CalculateElementForces ( ) const

This function calculates the internal element forces.

◆ CalculateGlobalNodalForces()

Vector Kratos::CrBeamElement3D2N::CalculateGlobalNodalForces ( ) const

◆ CalculateInitialLocalCS()

BoundedMatrix< double, CrBeamElement3D2N::msElementSize, CrBeamElement3D2N::msElementSize > Kratos::CrBeamElement3D2N::CalculateInitialLocalCS ( ) const

This function calculates the initial transformation matrix to globalize/localize vectors and/or matrices.

◆ CalculateLeftHandSide()

void Kratos::CrBeamElement3D2N::CalculateLeftHandSide ( MatrixType rLeftHandSideMatrix,
const ProcessInfo rCurrentProcessInfo 
)
overridevirtual

this is called during the assembling process in order to calculate the elemental left hand side matrix only

Parameters
rLeftHandSideMatrixthe elemental left hand side matrix
rCurrentProcessInfothe current process info instance

Reimplemented from Kratos::Element.

Reimplemented in Kratos::CrBeamElementLinear3D2N.

◆ CalculateLocalNodalForces()

Vector Kratos::CrBeamElement3D2N::CalculateLocalNodalForces ( ) const

This function calculates the local nodal forces.

Parameters
BisectrixThe bisectrix between the local axis1 from the last iter. step and the updated axis 1
VectorDifferenceThe vector differences of the quaternions

◆ CalculateLocalSystem()

void Kratos::CrBeamElement3D2N::CalculateLocalSystem ( MatrixType rLeftHandSideMatrix,
VectorType rRightHandSideVector,
const ProcessInfo rCurrentProcessInfo 
)
overridevirtual

ELEMENTS inherited from this class have to implement next CalculateLocalSystem, CalculateLeftHandSide and CalculateRightHandSide methods they can be managed internally with a private method to do the same calculations only once: MANDATORY this is called during the assembling process in order to calculate all elemental contributions to the global system matrix and the right hand side

Parameters
rLeftHandSideMatrixthe elemental left hand side matrix
rRightHandSideVectorthe elemental right hand side
rCurrentProcessInfothe current process info instance

Reimplemented from Kratos::Element.

Reimplemented in Kratos::CrBeamElementLinear3D2N.

◆ CalculateLumpedMassMatrix()

void Kratos::CrBeamElement3D2N::CalculateLumpedMassMatrix ( MatrixType rMassMatrix,
const ProcessInfo rCurrentProcessInfo 
) const

This function calculates the lumped mass matrix.

Parameters
rMassMatrixThe current mass matrix
rCurrentProcessInfoThe current Process information

◆ CalculateMassMatrix()

void Kratos::CrBeamElement3D2N::CalculateMassMatrix ( MatrixType rMassMatrix,
const ProcessInfo rCurrentProcessInfo 
)
overridevirtual

ELEMENTS inherited from this class must implement this methods if they need to add dynamic element contributions CalculateMassMatrix, CalculateDampingMatrix and CalculateLumpedMassVector methods are: OPTIONAL this is called during the assembling process in order to calculate the elemental mass matrix

Parameters
rMassMatrixthe elemental mass matrix
rCurrentProcessInfothe current process info instance

Reimplemented from Kratos::Element.

Reimplemented in Kratos::CrBeamElementLinear3D2N.

◆ CalculateOnIntegrationPoints()

void Kratos::CrBeamElement3D2N::CalculateOnIntegrationPoints ( const Variable< array_1d< double, 3 > > &  rVariable,
std::vector< array_1d< double, 3 > > &  rOutput,
const ProcessInfo rCurrentProcessInfo 
)
override

◆ CalculatePsi()

double Kratos::CrBeamElement3D2N::CalculatePsi ( const double  I,
const double  A_eff 
) const

This function calculates reduction values in case of shear-deformable structures.

Parameters
IThe second moment of area
A_effThe shear-effective area

◆ CalculateRightHandSide()

void Kratos::CrBeamElement3D2N::CalculateRightHandSide ( VectorType rRightHandSideVector,
const ProcessInfo rCurrentProcessInfo 
)
overridevirtual

this is called during the assembling process in order to calculate the elemental right hand side vector only

Parameters
rRightHandSideVectorthe elemental right hand side vector
rCurrentProcessInfothe current process info instance

Reimplemented from Kratos::Element.

Reimplemented in Kratos::CrBeamElementLinear3D2N.

◆ CalculateShearModulus()

double Kratos::CrBeamElement3D2N::CalculateShearModulus ( ) const

This function calculates shear modulus from user input values.

◆ CalculateSymmetricDeformationMode()

Vector Kratos::CrBeamElement3D2N::CalculateSymmetricDeformationMode ( ) const

This function calculates the symmetric deformation modes.

◆ CalculateTransformationS()

BoundedMatrix< double, CrBeamElement3D2N::msElementSize, CrBeamElement3D2N::msLocalSize > Kratos::CrBeamElement3D2N::CalculateTransformationS ( ) const

This function calculates a transformation matrix from deformation modes to real deformations.

◆ Check()

int Kratos::CrBeamElement3D2N::Check ( const ProcessInfo rCurrentProcessInfo) const
overridevirtual

This method provides the place to perform checks on the completeness of the input and the compatibility with the problem options as well as the contitutive laws selected It is designed to be called only once (or anyway, not often) typically at the beginning of the calculations, so to verify that nothing is missing from the input or that no common error is found.

Parameters
rCurrentProcessInfothis method is: MANDATORY

Reimplemented from Kratos::Element.

◆ ConstCalculateLeftHandSide()

void Kratos::CrBeamElement3D2N::ConstCalculateLeftHandSide ( MatrixType rLeftHandSideMatrix,
const ProcessInfo rCurrentProcessInfo 
) const

◆ ConstCalculateLocalSystem()

void Kratos::CrBeamElement3D2N::ConstCalculateLocalSystem ( MatrixType rLeftHandSideMatrix,
VectorType rRightHandSideVector,
const ProcessInfo rCurrentProcessInfo 
) const

◆ ConstCalculateRightHandSide()

void Kratos::CrBeamElement3D2N::ConstCalculateRightHandSide ( VectorType rRightHandSideVector,
const ProcessInfo rCurrentProcessInfo 
) const
virtual

Reimplemented in Kratos::GeoCrBeamElement3D2N.

◆ Create() [1/2]

Element::Pointer Kratos::CrBeamElement3D2N::Create ( IndexType  NewId,
GeometryType::Pointer  pGeom,
PropertiesType::Pointer  pProperties 
) const
overridevirtual

Creates a new element.

Parameters
NewIdThe Id of the new created element
pGeomThe pointer to the geometry of the element
pPropertiesThe pointer to property
Returns
The pointer to the created element

Reimplemented from Kratos::Element.

Reimplemented in Kratos::CrBeamElementLinear3D2N.

◆ Create() [2/2]

Element::Pointer Kratos::CrBeamElement3D2N::Create ( IndexType  NewId,
NodesArrayType const &  ThisNodes,
PropertiesType::Pointer  pProperties 
) const
overridevirtual

Creates a new element.

Parameters
NewIdThe Id of the new created element
ThisNodesThe array containing nodes
pPropertiesThe pointer to property
Returns
The pointer to the created element

Reimplemented from Kratos::Element.

Reimplemented in Kratos::CrBeamElementLinear3D2N.

◆ CreateElementStiffnessMatrix_Geometry()

BoundedMatrix< double, CrBeamElement3D2N::msElementSize, CrBeamElement3D2N::msElementSize > Kratos::CrBeamElement3D2N::CreateElementStiffnessMatrix_Geometry ( ) const

This function calculates the geometric part of the total stiffness matrix.

◆ CreateElementStiffnessMatrix_Material()

BoundedMatrix< double, CrBeamElement3D2N::msElementSize, CrBeamElement3D2N::msElementSize > Kratos::CrBeamElement3D2N::CreateElementStiffnessMatrix_Material ( ) const

This function calculates the elastic part of the total stiffness matrix.

◆ EquationIdVector()

void Kratos::CrBeamElement3D2N::EquationIdVector ( EquationIdVectorType rResult,
const ProcessInfo rCurrentProcessInfo 
) const
overridevirtual

ELEMENTS inherited from this class have to implement next EquationIdVector and GetDofList methods: MANDATORY this determines the elemental equation ID vector for all elemental DOFs

Parameters
rResultthe elemental equation ID vector
rCurrentProcessInfothe current process info instance

Reimplemented from Kratos::Element.

◆ FinalizeNonLinearIteration()

void Kratos::CrBeamElement3D2N::FinalizeNonLinearIteration ( const ProcessInfo rCurrentProcessInfo)
overridevirtual

this is called for non-linear analysis at the end of the iteration process

Reimplemented from Kratos::Element.

◆ GetCurrentNodalPosition()

BoundedVector< double, CrBeamElement3D2N::msLocalSize > Kratos::CrBeamElement3D2N::GetCurrentNodalPosition ( ) const

This function calculates the current nodal position.

◆ GetDofList()

void Kratos::CrBeamElement3D2N::GetDofList ( DofsVectorType rElementalDofList,
const ProcessInfo rCurrentProcessInfo 
) const
overridevirtual

determines the elemental list of DOFs

Parameters
ElementalDofListthe list of DOFs
rCurrentProcessInfothe current process info instance

Reimplemented from Kratos::Element.

◆ GetFirstDerivativesVector()

void Kratos::CrBeamElement3D2N::GetFirstDerivativesVector ( Vector values,
int  Step = 0 
) const
overridevirtual

Getting method to obtain the time derivative of variable which defines the degrees of freedom

Reimplemented from Kratos::Element.

◆ GetIncrementDeformation()

Vector Kratos::CrBeamElement3D2N::GetIncrementDeformation ( ) const

◆ GetIntegrationMethod()

CrBeamElement3D2N::IntegrationMethod Kratos::CrBeamElement3D2N::GetIntegrationMethod ( ) const
overridevirtual

returns the used integration method. In the general case this is the default integration method of the used geometry. I an other integration method is used the method has to be overwritten within the element

Returns
default integration method of the used Geometry

Reimplemented from Kratos::Element.

◆ GetSecondDerivativesVector()

void Kratos::CrBeamElement3D2N::GetSecondDerivativesVector ( Vector values,
int  Step = 0 
) const
overridevirtual

Getting method to obtain the second time derivative of variable which defines the degrees of freedom

Reimplemented from Kratos::Element.

◆ GetSpecifications()

const Parameters Kratos::CrBeamElement3D2N::GetSpecifications ( ) const
overridevirtual

This method provides the specifications/requirements of the element.

This can be used to enhance solvers and analysis. The following is an example: { "time_integration" : [], // NOTE: Options are static, implicit, explicit "framework" : "eulerian", // NOTE: Options are eulerian, lagrangian, ALE "symmetric_lhs" : true, // NOTE: Options are true/false "positive_definite_lhs" : false, // NOTE: Options are true/false "output" : { // NOTE: Values compatible as output "gauss_point" : ["INTEGRATION_WEIGTH"], "nodal_historical" : ["DISPLACEMENT"], "nodal_non_historical" : [], "entity" : [] }, "required_variables" : ["DISPLACEMENT"], // NOTE: Fill with the required variables "required_dofs" : ["DISPLACEMENT_X", "DISPLACEMENT_Y"], // NOTE: Fill with the required dofs "flags_used" : ["BOUNDARY", "ACTIVE"], // NOTE: Fill with the flags used "compatible_geometries" : ["Triangle2D3"], // NOTE: Compatible geometries. Options are "Point2D", "Point3D", "Sphere3D1", "Line2D2", "Line2D3", "Line3D2", "Line3D3", "Triangle2D3", "Triangle2D6", "Triangle3D3", "Triangle3D6", "Quadrilateral2D4", "Quadrilateral2D8", "Quadrilateral2D9", "Quadrilateral3D4", "Quadrilateral3D8", "Quadrilateral3D9", "Tetrahedra3D4" , "Tetrahedra3D10" , "Prism3D6" , "Prism3D15" , "Hexahedra3D8" , "Hexahedra3D20" , "Hexahedra3D27" "element_integrates_in_time" : true, // NOTE: Options are true/false "compatible_constitutive_laws": { "type" : ["PlaneStress","PlaneStrain"], // NOTE: List of CL compatible types. Options are "PlaneStress", "PlaneStrain", "3D" "dimension" : ["2D", "2D"], // NOTE: List of dimensions. Options are "2D", "3D", "2DAxysimm" "strain_size" : [3,3] // NOTE: List of strain sizes }, "documentation" : "This is an element" // NOTE: The documentation of the entity }

Returns
specifications The required specifications/requirements

Reimplemented from Kratos::Element.

◆ GetTransformationMatrixGlobal()

BoundedMatrix< double, CrBeamElement3D2N::msElementSize, CrBeamElement3D2N::msElementSize > Kratos::CrBeamElement3D2N::GetTransformationMatrixGlobal ( ) const

◆ GetValuesVector()

void Kratos::CrBeamElement3D2N::GetValuesVector ( Vector values,
int  Step = 0 
) const
overridevirtual

ELEMENTS inherited from this class must implement this methods if they need the values of the time derivatives of any of the dof set by the element. If the derivatives do not exist can set to zero these methods are: MANDATORY ( when compatibility with dynamics is required ) Getting method to obtain the variable which defines the degrees of freedom

Reimplemented from Kratos::Element.

◆ InitializeNonLinearIteration()

void Kratos::CrBeamElement3D2N::InitializeNonLinearIteration ( const ProcessInfo rCurrentProcessInfo)
overridevirtual

this is called for non-linear analysis at the beginning of the iteration process

Reimplemented from Kratos::Element.

◆ KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION()

Kratos::CrBeamElement3D2N::KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION ( CrBeamElement3D2N  )

◆ SaveQuaternionParameters()

void Kratos::CrBeamElement3D2N::SaveQuaternionParameters ( )

◆ UpdateQuaternionParameters()

void Kratos::CrBeamElement3D2N::UpdateQuaternionParameters ( double rScalNodeA,
double rScalNodeB,
Vector rVecNodeA,
Vector rVecNodeB 
) const

◆ UpdateRotationMatrixLocal()

BoundedMatrix< double, CrBeamElement3D2N::msDimension, CrBeamElement3D2N::msDimension > Kratos::CrBeamElement3D2N::UpdateRotationMatrixLocal ( Vector Bisectrix,
Vector VectorDifference 
) const

This function updates constantly the transformation matrix.

Friends And Related Function Documentation

◆ Serializer

friend class Serializer
friend

Member Data Documentation

◆ msDimension

constexpr int Kratos::CrBeamElement3D2N::msDimension = 3
staticconstexprprotected

◆ msElementSize

constexpr unsigned int Kratos::CrBeamElement3D2N::msElementSize = msLocalSize * 2
staticconstexprprotected

◆ msLocalSize

constexpr unsigned int Kratos::CrBeamElement3D2N::msLocalSize = msNumberOfNodes * msDimension
staticconstexprprotected

◆ msNumberOfNodes

constexpr int Kratos::CrBeamElement3D2N::msNumberOfNodes = 2
staticconstexprprotected

The documentation for this class was generated from the following files: