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.
ring_element_3D.hpp
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Klaus B. Sautter
11 //
12 //
13 //
14 
15 #if !defined(KRATOS_RING_ELEMENT_3D_H_INCLUDED )
16 #define KRATOS_RING_ELEMENT_3D_H_INCLUDED
17 
18 // System includes
19 
20 // External includes
21 
22 // Project includes
23 #include "includes/element.h"
24 #include "includes/define.h"
25 
26 namespace Kratos
27 {
36  class KRATOS_API(CABLE_NET_APPLICATION) RingElement3D : public Element
37  {
38  protected:
39 
40  public:
42 
43 
44  typedef Element BaseType;
52  typedef BaseType::EquationIdVectorType EquationIdVectorType;
53  typedef BaseType::DofsVectorType DofsVectorType;
54 
55 
58  GeometryType::Pointer pGeometry);
60  GeometryType::Pointer pGeometry,
61  PropertiesType::Pointer pProperties);
62 
63 
64  ~RingElement3D() override;
65 
73  Element::Pointer Create(
74  IndexType NewId,
75  GeometryType::Pointer pGeom,
76  PropertiesType::Pointer pProperties
77  ) const override;
78 
86  Element::Pointer Create(
87  IndexType NewId,
88  NodesArrayType const& ThisNodes,
89  PropertiesType::Pointer pProperties
90  ) const override;
91 
92  void EquationIdVector(
93  EquationIdVectorType& rResult,
94  const ProcessInfo& rCurrentProcessInfo) const override;
95 
96  void GetDofList(
97  DofsVectorType& rElementalDofList,
98  const ProcessInfo& rCurrentProcessInfo) const override;
99 
100  void CalculateLeftHandSide(
101  MatrixType& rLeftHandSideMatrix,
102  const ProcessInfo& rCurrentProcessInfo) override;
103 
104  void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix,
105  VectorType &rRightHandSideVector,
106  const ProcessInfo &rCurrentProcessInfo) override;
107 
108  void CalculateRightHandSide(VectorType &rRightHandSideVector,
109  const ProcessInfo &rCurrentProcessInfo) override;
110 
111  void GetValuesVector(Vector& rValues,int Step = 0) const override;
112  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
113  void GetFirstDerivativesVector(Vector& rValues,int Step = 0) const override;
114 
115  Matrix ElasticStiffnessMatrix() const;
116  Matrix GeometricStiffnessMatrix() const;
117  inline Matrix TotalStiffnessMatrix() const;
118 
119  double GetCurrentLength() const;
120  double GetRefLength() const;
121  double CalculateGreenLagrangeStrain() const;
122  double LinearStiffness() const;
123 
124  void CalculateLumpedMassVector(
125  VectorType &rLumpedMassVector,
126  const ProcessInfo& rCurrentProcessInfo) const override;
127 
128  void CalculateMassMatrix(
129  MatrixType& rMassMatrix,
130  const ProcessInfo& rCurrentProcessInfo) override;
131 
132  void CalculateDampingMatrix(
133  MatrixType& rDampingMatrix,
134  const ProcessInfo& rCurrentProcessInfo) override;
135 
136 
137  void AddExplicitContribution(
138  const VectorType& rRHSVector,
139  const Variable<VectorType>& rRHSVariable,
140  const Variable<double >& rDestinationVariable,
141  const ProcessInfo& rCurrentProcessInfo
142  ) override;
143 
144  void AddExplicitContribution(const VectorType& rRHSVector,
145  const Variable<VectorType>& rRHSVariable,
146  const Variable<array_1d<double, 3> >& rDestinationVariable,
147  const ProcessInfo& rCurrentProcessInfo
148  ) override;
149 
150  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
151 
155  bool HasSelfWeight() const;
156 
160  Vector CalculateBodyForces();
161 
162  private:
163 
164  Vector GetCurrentLengthArray() const;
165  Vector GetRefLengthArray() const;
166  Vector GetDeltaPositions(const int& rDirection) const;
167  Vector GetDirectionVectorNt() const;
168  Vector GetInternalForces() const;
169 
170  friend class Serializer;
171  void save(Serializer& rSerializer) const override;
172  void load(Serializer& rSerializer) override;
173  };
174 
175 
176 }
177 
178 
179 #endif
Base class for all Elements.
Definition: element.h:60
std::size_t IndexType
Defines the index type.
Definition: geometrical_object.h:73
Geometry base class.
Definition: geometry.h:71
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
This is a ring elemen with 3 translational dofs per node.
Definition: ring_element_3D.hpp:37
Element BaseType
Definition: ring_element_3D.hpp:44
BaseType::MatrixType MatrixType
Definition: ring_element_3D.hpp:50
BaseType::EquationIdVectorType EquationIdVectorType
Definition: ring_element_3D.hpp:52
RingElement3D()
Definition: ring_element_3D.hpp:56
BaseType::PropertiesType PropertiesType
Definition: ring_element_3D.hpp:47
BaseType::IndexType IndexType
Definition: ring_element_3D.hpp:48
BaseType::DofsVectorType DofsVectorType
Definition: ring_element_3D.hpp:53
BaseType::NodesArrayType NodesArrayType
Definition: ring_element_3D.hpp:46
BaseType::GeometryType GeometryType
Definition: ring_element_3D.hpp:45
BaseType::SizeType SizeType
Definition: ring_element_3D.hpp:49
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(RingElement3D)
BaseType::VectorType VectorType
Definition: ring_element_3D.hpp:51
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 rRH...
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
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
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