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.
cr_beam_element_2D2N.hpp
Go to the documentation of this file.
1 // KRATOS ___| | | |
2 // \___ \ __| __| | | __| __| | | __| _` | |
3 // | | | | | ( | | | | ( | |
4 // _____/ \__|_| \__,_|\___|\__|\__,_|_| \__,_|_| MECHANICS
5 //
6 // License: BSD License
7 // license: StructuralMechanicsApplication/license.txt
8 //
9 // Main authors: Klaus B. Sautter
10 //
11 
12 #pragma once
13 
14 // System includes
15 
16 // External includes
17 
18 // Project includes
19 #include "includes/element.h"
20 #include "includes/define.h"
21 #include "includes/variables.h"
22 #include "includes/serializer.h"
23 
24 namespace Kratos
25 {
34 class KRATOS_API(STRUCTURAL_MECHANICS_APPLICATION) CrBeamElement2D2N : public Element
35 {
36 protected:
37  //const values
38  static constexpr int msNumberOfNodes = 2;
39  static constexpr int msDimension = 2;
40  static constexpr unsigned int msLocalSize = 3;
41  static constexpr unsigned int msElementSize = msLocalSize * 2;
42 
43  // stores the deformation modes
44  BoundedVector<double,msLocalSize> mDeformationForces = ZeroVector(msLocalSize);
45 
46 
47  // stores the globalized internal forces for calculation of the residual
48  Vector mInternalGlobalForces = ZeroVector(msElementSize);
49 
50 public:
52 
53 
54  typedef Element BaseType;
62  typedef BaseType::EquationIdVectorType EquationIdVectorType;
63  typedef BaseType::DofsVectorType DofsVectorType;
64 
66  CrBeamElement2D2N(IndexType NewId, GeometryType::Pointer pGeometry);
67  CrBeamElement2D2N(IndexType NewId, GeometryType::Pointer pGeometry,
68  PropertiesType::Pointer pProperties);
69 
70 
71  ~CrBeamElement2D2N() override;
72 
73 
81  Element::Pointer Create(
82  IndexType NewId,
83  GeometryType::Pointer pGeom,
84  PropertiesType::Pointer pProperties
85  ) const override;
86 
94  Element::Pointer Create(
95  IndexType NewId,
96  NodesArrayType const& ThisNodes,
97  PropertiesType::Pointer pProperties
98  ) const override;
99 
100  void EquationIdVector(
101  EquationIdVectorType& rResult,
102  const ProcessInfo& rCurrentProcessInfo) const override;
103 
104  void GetDofList(
105  DofsVectorType& rElementalDofList,
106  const ProcessInfo& rCurrentProcessInfo) const override;
107 
108  void GetValuesVector(
109  Vector& rValues,
110  int Step = 0) const override;
111 
112  void GetSecondDerivativesVector(
113  Vector& rValues,
114  int Step = 0) const override;
115 
116  void GetFirstDerivativesVector(
117  Vector& rValues,
118  int Step = 0) const override;
119 
120  void CalculateMassMatrix(
121  MatrixType& rMassMatrix,
122  const ProcessInfo& rCurrentProcessInfo) override;
123 
124  void CalculateDampingMatrix(
125  MatrixType& rDampingMatrix,
126  const ProcessInfo& rCurrentProcessInfo) override;
127 
128  void CalculateLocalSystem(
129  MatrixType& rLeftHandSideMatrix,
130  VectorType& rRightHandSideVector,
131  const ProcessInfo& rCurrentProcessInfo) override;
132 
133  void CalculateRightHandSide(
134  VectorType& rRightHandSideVector,
135  const ProcessInfo& rCurrentProcessInfo) override;
136 
137  void CalculateLeftHandSide(
138  MatrixType& rLeftHandSideMatrix,
139  const ProcessInfo& rCurrentProcessInfo) override;
140 
141  void AddExplicitContribution(const VectorType& rRHSVector,
142  const Variable<VectorType>& rRHSVariable,
143  const Variable<array_1d<double, 3> >& rDestinationVariable,
144  const ProcessInfo& rCurrentProcessInfo) override;
145 
146  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
147 
151 
155  double CalculateShearModulus() const;
156 
157 
163  double CalculatePsi(const double I, const double A_eff) const;
164 
168  double CalculateInitialElementAngle() const;
169 
173  double CalculateDeformedElementAngle();
174 
178  BoundedVector<double,msElementSize> CalculateBodyForces() const;
179 
186  void CalculateAndAddWorkEquivalentNodalForcesLineLoad(
187  const BoundedVector<double,3> ForceInput,
188  BoundedVector<double,msElementSize>& rRightHandSideVector,
189  const double GeometryLength)const ;
190 
191  IntegrationMethod GetIntegrationMethod() const override;
192 
196  BoundedMatrix<double,msElementSize,msLocalSize> CalculateTransformationS() const;
197 
201  virtual double CalculateLength() const;
202 
206  BoundedMatrix<double,msLocalSize,msLocalSize> CreateElementStiffnessMatrix_Kd_mat() const;
207 
211  BoundedMatrix<double,msLocalSize,msLocalSize> CreateElementStiffnessMatrix_Kd_geo() const;
212 
216  BoundedMatrix<double,msElementSize,msElementSize> CreateElementStiffnessMatrix_Kr() const;
217 
221  BoundedMatrix<double,msElementSize,msElementSize> CreateElementStiffnessMatrix_Total() const;
222 
223 
228  void GlobalizeMatrix(Matrix& A);
229 
234  void GlobalizeVector(Vector& A);
235 
240  double Modulus2Pi(double A) const;
241 
245  virtual BoundedMatrix<double,msElementSize,msElementSize> CreateRotationMatrix();
246 
247 
251  BoundedVector<double,msLocalSize> CalculateDeformationParameters();
252 
256  BoundedVector<double,msLocalSize> CalculateInternalStresses_DeformationModes();
257 
261  BoundedVector<double,msElementSize> ReturnElementForces_Local();
262 
263 
269  const Variable<array_1d<double, 3 > >& rVariable,
270  std::vector< array_1d<double, 3 > >& rOutput,
271  const ProcessInfo& rCurrentProcessInfo) override;
272 
273  const Parameters GetSpecifications() const override;
274 
275 private:
276 
277  friend class Serializer;
278  void save(Serializer& rSerializer) const override;
279  void load(Serializer& rSerializer) override;
280 
281 };
282 
283 }
This is a 2D-2node beam element with 2 translational dofs and 1 rotational dof per node.
Definition: cr_beam_element_2D2N.hpp:35
BaseType::SizeType SizeType
Definition: cr_beam_element_2D2N.hpp:59
BaseType::DofsVectorType DofsVectorType
Definition: cr_beam_element_2D2N.hpp:63
BaseType::NodesArrayType NodesArrayType
Definition: cr_beam_element_2D2N.hpp:56
Element BaseType
Definition: cr_beam_element_2D2N.hpp:54
BaseType::EquationIdVectorType EquationIdVectorType
Definition: cr_beam_element_2D2N.hpp:62
BaseType::PropertiesType PropertiesType
Definition: cr_beam_element_2D2N.hpp:57
BaseType::GeometryType GeometryType
Definition: cr_beam_element_2D2N.hpp:55
CrBeamElement2D2N()
Definition: cr_beam_element_2D2N.hpp:65
BaseType::VectorType VectorType
Definition: cr_beam_element_2D2N.hpp:61
BaseType::IndexType IndexType
Definition: cr_beam_element_2D2N.hpp:58
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(CrBeamElement2D2N)
BaseType::MatrixType MatrixType
Definition: cr_beam_element_2D2N.hpp:60
Base class for all Elements.
Definition: element.h:60
std::size_t IndexType
Defines the index type.
Definition: geometrical_object.h:73
IntegrationMethod
Definition: geometry_data.h:76
Geometry base class.
Definition: geometry.h:71
Definition: amatrix_interface.h:41
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
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
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
Vector VectorType
Definition: geometrical_transformation_utilities.h:56
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
std::size_t SizeType
Definition: nurbs_utilities.h:41
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
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
Properties PropertiesType
Definition: regenerate_pfem_pressure_conditions_process.h:26
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
def load(f)
Definition: ode_solve.py:307
A
Definition: sensitivityMatrix.py:70