13 #if !defined(KRATOS_GEO_CURVED_BEAM_ELEMENT_H_INCLUDED)
14 #define KRATOS_GEO_CURVED_BEAM_ELEMENT_H_INCLUDED
37 template <
unsigned int TDim,
unsigned int TNumNodes>
85 mThisIntegrationMethod = this->GetIntegrationMethod();
95 PropertiesType::Pointer pProperties)
const override;
97 Element::Pointer
Create(
IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties)
const override;
99 int Check(
const ProcessInfo& rCurrentProcessInfo)
const override;
106 std::vector<Matrix>& rOutput,
118 SizeType GetCrossNumberIntegrationPoints()
const override;
119 SizeType GetAlongNumberIntegrationPoints()
const override;
121 double CalculateIntegrationCoefficient(
unsigned int GPointCross,
double detJ,
double weight)
const;
123 virtual void CalculateBMatrix(
Matrix&
B,
124 unsigned int GPointCross,
126 ElementVariables& rVariables)
const;
128 virtual void CalculateLocalBMatrix(
Matrix&
B,
129 unsigned int GPointCross,
131 ElementVariables& rVariables)
const;
133 void CalculateStrainVector(ElementVariables& rVariables)
const;
136 virtual void CalculateAll(
MatrixType& rLeftHandSideMatrix,
139 const bool CalculateStiffnessMatrixFlag,
140 const bool CalculateResidualVectorFlag)
override;
142 virtual void CalculateAndAddLHS(
MatrixType& rLeftHandSideMatrix, ElementVariables& rVariables)
const;
144 virtual void CalculateAndAddRHS(
VectorType& rRightHandSideVector,
145 ElementVariables& rVariables,
146 unsigned int GPoint)
const;
148 virtual void CalculateLocalInternalForce(
VectorType& rRightHandSideVector,
const ProcessInfo& rCurrentProcessInfo);
152 virtual void CalculateNodalCrossDirection(
Matrix& NodalCrossDirection)
const override;
154 virtual double CalculateAngleAtGaussPoint(
const Matrix& GradNe)
const;
156 virtual double CalculateAngleAtNode(
unsigned int GPoint,
159 virtual void CalculateJacobianMatrix(
unsigned int GPointCross,
165 void CalculateAndAddStiffnessForce(
VectorType& rRightHandSideVector,
167 unsigned int GPoint)
const;
174 const ProcessInfo& rCurrentProcessInfo)
const override;
176 void InterpolateOnOutputPoints(
Matrix& Values)
const;
177 void InterpolateOnOutputPoints(
Vector& Values)
const;
192 void save(
Serializer& rSerializer)
const override
203 static constexpr
SizeType N_DOF_NODE_DISP = TDim;
204 static constexpr
SizeType N_DOF_NODE_ROT = (TDim == 2 ? 1 : 3);
205 static constexpr
SizeType N_DOF_NODE = N_DOF_NODE_DISP + N_DOF_NODE_ROT;
206 static constexpr
SizeType N_POINT_CROSS = 2;
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Base class for all Elements.
Definition: element.h:60
std::size_t SizeType
Definition: element.h:94
std::size_t IndexType
Definition: flags.h:74
This is a geometrically non-linear (curved) beam element. The formulation can be found in papers writ...
Definition: geo_curved_beam_element.hpp:40
KRATOS_CLASS_POINTER_DEFINITION(GeoCurvedBeamElement)
GeoCurvedBeamElement(IndexType NewId=0)
Default Constructor.
Definition: geo_curved_beam_element.hpp:67
GeoCurvedBeamElement(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor using Geometry.
Definition: geo_curved_beam_element.hpp:76
typename GeoStructuralBaseElement< TDim, TNumNodes >::ElementVariables ElementVariables
Definition: geo_curved_beam_element.hpp:60
GeoCurvedBeamElement(IndexType NewId, const NodesArrayType &ThisNodes)
Constructor using an array of nodes.
Definition: geo_curved_beam_element.hpp:70
void SetRotationalInertiaVector(const PropertiesType &Prop, Vector &RotationalInertia) const
virtual ~GeoCurvedBeamElement()
Destructor.
Definition: geo_curved_beam_element.hpp:89
GeoCurvedBeamElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Constructor using Properties.
Definition: geo_curved_beam_element.hpp:82
virtual void CalculateTransformationMatrix(Matrix &TransformationMatrix, const Matrix &GradNe) const
Definition: geo_structural_base_element.hpp:28
DenseVector< Matrix > ShapeFunctionsGradientsType
Definition: geometry_data.h:208
Geometry base class.
Definition: geometry.h:71
PointerVector< TPointType > PointsArrayType
Definition: geometry.h:118
This class defines the node.
Definition: node.h:65
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
pybind11::list CalculateOnIntegrationPoints(TObject &dummy, const Variable< TDataType > &rVariable, const ProcessInfo &rProcessInfo)
Definition: add_mesh_to_python.cpp:142
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
Internals::Matrix< double, AMatrix::dynamic, 1 > Vector
Definition: amatrix_interface.h:472
Internals::Matrix< double, AMatrix::dynamic, AMatrix::dynamic > Matrix
Definition: amatrix_interface.h:470
std::size_t SizeType
The definition of the size type.
Definition: mortar_classes.h:43
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
def load(f)
Definition: ode_solve.py:307
B
Definition: sensitivityMatrix.py:76
Definition: constitutive_law.h:189
Member Variables.
Definition: geo_structural_base_element.hpp:114