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.
empirical_spring.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_EMPIRICAL_SPRING_H_INCLUDED )
16 #define KRATOS_EMPIRICAL_SPRING_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 {
37  class KRATOS_API(CABLE_NET_APPLICATION) EmpiricalSpringElement3D2N : public Element
38  {
39  protected:
40  //const values
41  static constexpr int msNumberOfNodes = 2;
42  static constexpr int msDimension = 3;
43  static constexpr unsigned int msLocalSize = msNumberOfNodes * msDimension;
44 
45  public:
47 
48 
49  typedef Element BaseType;
57  typedef BaseType::EquationIdVectorType EquationIdVectorType;
58  typedef BaseType::DofsVectorType DofsVectorType;
59 
60 
63  GeometryType::Pointer pGeometry);
65  GeometryType::Pointer pGeometry,
66  PropertiesType::Pointer pProperties);
67 
68 
69  ~EmpiricalSpringElement3D2N() override;
70 
78  Element::Pointer Create(
79  IndexType NewId,
80  GeometryType::Pointer pGeom,
81  PropertiesType::Pointer pProperties
82  ) const override;
83 
91  Element::Pointer Create(
92  IndexType NewId,
93  NodesArrayType const& ThisNodes,
94  PropertiesType::Pointer pProperties
95  ) const override;
96 
97  void EquationIdVector(
98  EquationIdVectorType& rResult,
99  const ProcessInfo& rCurrentProcessInfo) const override;
100 
101  void GetDofList(
102  DofsVectorType& rElementalDofList,
103  const ProcessInfo& rCurrentProcessInfo) const override;
104 
109  CreateElementStiffnessMatrix(const ProcessInfo& rCurrentProcessInfo);
110 
111 
112  void CalculateLocalSystem(
113  MatrixType& rLeftHandSideMatrix,
114  VectorType& rRightHandSideVector,
115  const ProcessInfo& rCurrentProcessInfo) override;
116 
117 
118  void CalculateRightHandSide(
119  VectorType& rRightHandSideVector,
120  const ProcessInfo& rCurrentProcessInfo) override;
121 
122  void CalculateLeftHandSide(
123  MatrixType& rLeftHandSideMatrix,
124  const ProcessInfo& rCurrentProcessInfo) override;
125 
126  void AddExplicitContribution(
127  const VectorType& rRHSVector, const Variable<VectorType>& rRHSVariable,
128  const Variable<array_1d<double, 3>>& rDestinationVariable,
129  const ProcessInfo& rCurrentProcessInfo) override;
130 
131  void AddExplicitContribution(
132  const VectorType& rRHSVector,
133  const Variable<VectorType>& rRHSVariable,
134  const Variable<double >& rDestinationVariable,
135  const ProcessInfo& rCurrentProcessInfo) override;
136 
137  void CalculateMassMatrix(
138  MatrixType& rMassMatrix,
139  const ProcessInfo& rCurrentProcessInfo) override;
140 
141  void CalculateLumpedMassVector(
142  VectorType &rLumpedMassVector,
143  const ProcessInfo &rCurrentProcessInfo) const override;
144 
145  void CalculateDampingMatrix(MatrixType& rDampingMatrix,
146  const ProcessInfo& rCurrentProcessInfo) override;
147 
152  double EvaluatePolynomial(const Vector& rPolynomial) const;
153  double EvaluatePolynomialFirstDerivative(const Vector& rPolynomial) const;
154 
155  void LocalizeVector(BoundedVector<double,msLocalSize>& rInputVector);
156  void GlobalizeVector(BoundedVector<double,msLocalSize>& rInputVector);
157  void GlobalizeMatrix(BoundedMatrix<double,msLocalSize,msLocalSize>& rInputMatrix);
158 
159  void WriteTransformationCoordinates(BoundedVector<double,msLocalSize>
160  & rReferenceCoordinates);
161 
162  void CreateTransformationMatrix(BoundedMatrix<double,msLocalSize,
163  msLocalSize>& rRotationMatrix);
164 
165  double GetElementElongation() const;
166 
167 
168  void GetValuesVector(
169  Vector& rValues,
170  int Step = 0) const override;
171 
172  void GetSecondDerivativesVector(
173  Vector& rValues,
174  int Step = 0) const override;
175 
176  void GetFirstDerivativesVector(
177  Vector& rValues,
178  int Step = 0) const override;
179 
180  int Check(
181  const ProcessInfo& rCurrentProcessInfo) const override;
182 
183 private:
184 
185  friend class Serializer;
186  void save(Serializer& rSerializer) const override;
187  void load(Serializer& rSerializer) override;
188  };
189 
190 
191 }
192 
193 
194 #endif
Base class for all Elements.
Definition: element.h:60
This spring reads a fitted polynomial as displacement-load curve takes u and returns f.
Definition: empirical_spring.hpp:38
BaseType::PropertiesType PropertiesType
Definition: empirical_spring.hpp:52
BaseType::VectorType VectorType
Definition: empirical_spring.hpp:56
BaseType::DofsVectorType DofsVectorType
Definition: empirical_spring.hpp:58
BaseType::MatrixType MatrixType
Definition: empirical_spring.hpp:55
BaseType::IndexType IndexType
Definition: empirical_spring.hpp:53
Element BaseType
Definition: empirical_spring.hpp:49
BaseType::NodesArrayType NodesArrayType
Definition: empirical_spring.hpp:51
BaseType::GeometryType GeometryType
Definition: empirical_spring.hpp:50
EmpiricalSpringElement3D2N()
Definition: empirical_spring.hpp:61
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(EmpiricalSpringElement3D2N)
BaseType::SizeType SizeType
Definition: empirical_spring.hpp:54
BaseType::EquationIdVectorType EquationIdVectorType
Definition: empirical_spring.hpp:57
std::size_t IndexType
Defines the index type.
Definition: geometrical_object.h:73
Geometry base class.
Definition: geometry.h:71
Definition: amatrix_interface.h:41
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
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