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.
geo_curved_beam_element.hpp
Go to the documentation of this file.
1 // KRATOS___
2 // // ) )
3 // // ___ ___
4 // // ____ //___) ) // ) )
5 // // / / // // / /
6 // ((____/ / ((____ ((___/ / MECHANICS
7 //
8 // License: geo_mechanics_application/license.txt
9 //
10 // Main authors: Vahid Galavi
11 //
12 
13 #if !defined(KRATOS_GEO_CURVED_BEAM_ELEMENT_H_INCLUDED)
14 #define KRATOS_GEO_CURVED_BEAM_ELEMENT_H_INCLUDED
15 
16 // Project includes
17 
18 // Application includes
20 
21 namespace Kratos
22 {
37 template <unsigned int TDim, unsigned int TNumNodes>
38 class KRATOS_API(GEO_MECHANICS_APPLICATION) GeoCurvedBeamElement
39  : public GeoStructuralBaseElement<TDim, TNumNodes>
40 {
41 public:
42  using IndexType = std::size_t;
44  using NodeType = Node;
47  using VectorType = Vector;
48  using MatrixType = Matrix;
50 
52  using SizeType = std::size_t;
53 
61 
63 
65 
67  GeoCurvedBeamElement(IndexType NewId = 0) : GeoStructuralBaseElement<TDim, TNumNodes>(NewId) {}
68 
71  : GeoStructuralBaseElement<TDim, TNumNodes>(NewId, ThisNodes)
72  {
73  }
74 
76  GeoCurvedBeamElement(IndexType NewId, GeometryType::Pointer pGeometry)
77  : GeoStructuralBaseElement<TDim, TNumNodes>(NewId, pGeometry)
78  {
79  }
80 
82  GeoCurvedBeamElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
83  : GeoStructuralBaseElement<TDim, TNumNodes>(NewId, pGeometry, pProperties)
84  {
85  mThisIntegrationMethod = this->GetIntegrationMethod();
86  }
87 
89  virtual ~GeoCurvedBeamElement() {}
90 
92 
93  Element::Pointer Create(IndexType NewId,
94  NodesArrayType const& ThisNodes,
95  PropertiesType::Pointer pProperties) const override;
96 
97  Element::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override;
98 
99  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
100 
102 
103  void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) override;
104 
105  void CalculateOnIntegrationPoints(const Variable<Matrix>& rVariable,
106  std::vector<Matrix>& rOutput,
107  const ProcessInfo& rCurrentProcessInfo) override;
108 
110  std::vector<array_1d<double, 3>>& rOutput,
111  const ProcessInfo& rCurrentProcessInfo) override;
112 
114 
115 protected:
117 
118  SizeType GetCrossNumberIntegrationPoints() const override;
119  SizeType GetAlongNumberIntegrationPoints() const override;
120 
121  double CalculateIntegrationCoefficient(unsigned int GPointCross, double detJ, double weight) const;
122 
123  virtual void CalculateBMatrix(Matrix& B,
124  unsigned int GPointCross,
125  const BoundedMatrix<double, TDim, TDim>& InvertDetJacobian,
126  ElementVariables& rVariables) const;
127 
128  virtual void CalculateLocalBMatrix(Matrix& B,
129  unsigned int GPointCross,
130  const BoundedMatrix<double, TDim, TDim>& InvertDetJacobian,
131  ElementVariables& rVariables) const;
132 
133  void CalculateStrainVector(ElementVariables& rVariables) const;
134 
136  virtual void CalculateAll(MatrixType& rLeftHandSideMatrix,
137  VectorType& rRightHandSideVector,
138  const ProcessInfo& rCurrentProcessInfo,
139  const bool CalculateStiffnessMatrixFlag,
140  const bool CalculateResidualVectorFlag) override;
141 
142  virtual void CalculateAndAddLHS(MatrixType& rLeftHandSideMatrix, ElementVariables& rVariables) const;
143 
144  virtual void CalculateAndAddRHS(VectorType& rRightHandSideVector,
145  ElementVariables& rVariables,
146  unsigned int GPoint) const;
147 
148  virtual void CalculateLocalInternalForce(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo);
149 
150  virtual void CalculateTransformationMatrix(Matrix& TransformationMatrix, const Matrix& GradNe) const;
151 
152  virtual void CalculateNodalCrossDirection(Matrix& NodalCrossDirection) const override;
153 
154  virtual double CalculateAngleAtGaussPoint(const Matrix& GradNe) const;
155 
156  virtual double CalculateAngleAtNode(unsigned int GPoint,
157  const BoundedMatrix<double, TNumNodes, TNumNodes>& DN_DXContainer) const;
158 
159  virtual void CalculateJacobianMatrix(unsigned int GPointCross,
160  const ElementVariables& rVariables,
161  BoundedMatrix<double, TDim, TDim>& DeterminantJacobian) const;
162 
163  void CalculateAndAddBodyForce(VectorType& rRightHandSideVector, ElementVariables& rVariables) const;
164 
165  void CalculateAndAddStiffnessForce(VectorType& rRightHandSideVector,
166  ElementVariables& rVariables,
167  unsigned int GPoint) const;
168  void SetRotationalInertiaVector(const PropertiesType& Prop, Vector& RotationalInertia) const;
169 
170  void InitializeElementVariables(ElementVariables& rVariables,
171  ConstitutiveLaw::Parameters& rConstitutiveParameters,
172  const GeometryType& Geom,
173  const PropertiesType& Prop,
174  const ProcessInfo& rCurrentProcessInfo) const override;
175 
176  void InterpolateOnOutputPoints(Matrix& Values) const;
177  void InterpolateOnOutputPoints(Vector& Values) const;
178 
180 
181 private:
184 
187 
189 
190  friend class Serializer;
191 
192  void save(Serializer& rSerializer) const override
193  {
195  }
196 
197  void load(Serializer& rSerializer) override
198  {
200  }
201 
202  // const values
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;
207 
208 }; // Class GeoCurvedBeamElement
209 
210 } // namespace Kratos
211 
212 #endif // KRATOS_GEO_CURVED_BEAM_ELEMENT_H_INCLUDED defined
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
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
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