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_truss_element_base.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: Klaus B. Sautter,
11 // Vahid Galavi
12 //
13 
14 #if !defined(KRATOS_GEO_TRUSS_ELEMENT_BASE_H_INCLUDED)
15 #define KRATOS_GEO_TRUSS_ELEMENT_BASE_H_INCLUDED
16 
17 // System includes
18 
19 // External includes
20 
21 // Project includes
22 #include "includes/define.h"
23 #include "includes/element.h"
24 #include "includes/variables.h"
25 
26 namespace Kratos
27 {
35 template <unsigned int TDim, unsigned int TNumNodes>
36 class KRATOS_API(GEO_MECHANICS_APPLICATION) GeoTrussElementBase : public Element
37 {
38 protected:
39  ConstitutiveLaw::Pointer mpConstitutiveLaw = nullptr;
40  static constexpr SizeType NDof = TDim * TNumNodes;
41  static constexpr SizeType DIM = TDim;
42  static constexpr SizeType NUM_NODES = TNumNodes;
43 
44 public:
46 
47  using BaseType = Element;
55  using EquationIdVectorType = BaseType::EquationIdVectorType;
56  using DofsVectorType = BaseType::DofsVectorType;
59 
61  GeoTrussElementBase(IndexType NewId, GeometryType::Pointer pGeometry);
62  GeoTrussElementBase(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
63 
64  ~GeoTrussElementBase() override;
65 
73  Element::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override;
74 
82  Element::Pointer Create(IndexType NewId,
83  NodesArrayType const& ThisNodes,
84  PropertiesType::Pointer pProperties) const override;
85 
86  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
87 
88  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
89 
90  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
91 
95  virtual void CreateElementStiffnessMatrix(MatrixType& rLocalStiffnessMatrix,
96  const ProcessInfo& rCurrentProcessInfo);
97 
98  void Calculate(const Variable<Matrix>& rVariable, Matrix& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
99 
100  void Calculate(const Variable<double>& rVariable, double& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
101 
102  void CalculateOnIntegrationPoints(const Variable<double>& rVariable,
103  std::vector<double>& rOutput,
104  const ProcessInfo& rCurrentProcessInfo) override;
105 
107  std::vector<array_1d<double, 3>>& rOutput,
108  const ProcessInfo& rCurrentProcessInfo) override;
109 
110  void CalculateOnIntegrationPoints(const Variable<Vector>& rVariable,
111  std::vector<Vector>& rOutput,
112  const ProcessInfo& rCurrentProcessInfo) override;
113 
118  virtual void UpdateInternalForces(FullDofVectorType& rInternalForces, const ProcessInfo& rCurrentProcessInfo);
119 
125 
126  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
127  VectorType& rRightHandSideVector,
128  const ProcessInfo& rCurrentProcessInfo) override;
129 
130  void CalculateRightHandSide(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
131 
132  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo) override;
133 
134  void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) override;
135 
136  void CalculateConsistentMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) const;
137 
138  void CalculateDampingMatrix(MatrixType& rDampingMatrix, const ProcessInfo& rCurrentProcessInfo) override;
139 
149  void AddExplicitContribution(const VectorType& rRHSVector,
150  const Variable<VectorType>& rRHSVariable,
151  const Variable<double>& rDestinationVariable,
152  const ProcessInfo& rCurrentProcessInfo) override;
153 
163  void AddExplicitContribution(const VectorType& rRHSVector,
164  const Variable<VectorType>& rRHSVariable,
165  const Variable<array_1d<double, 3>>& rDestinationVariable,
166  const ProcessInfo& rCurrentProcessInfo) override;
167 
168  void GetValuesVector(Vector& rValues, int Step = 0) const override;
169 
170  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
171 
172  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
173 
174  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
175 
179  double CalculateGreenLagrangeStrain() const;
180 
184  void CalculateBodyForces(FullDofVectorType& rGlobalBodyForces);
185 
192  const ProcessInfo& rCurrentProcessInfo);
193 
199  void CalculateElasticStiffnessMatrix(MatrixType& rElasticStiffnessMatrix, const ProcessInfo& rCurrentProcessInfo);
200 
205  virtual void WriteTransformationCoordinates(FullDofVectorType& rReferenceCoordinates);
206 
207  double ReturnTangentModulus1D(const ProcessInfo& rCurrentProcessInfo);
208 
209  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
210 
211 private:
216  void CalculateLumpedMassVector(VectorType& rMassVector, const ProcessInfo& rCurrentProcessInfo) const override;
217 
218  friend class Serializer;
219  void save(Serializer& rSerializer) const override;
220  void load(Serializer& rSerializer) override;
221 };
222 
223 } // namespace Kratos
224 
225 #endif
Base class for all Elements.
Definition: element.h:60
std::size_t SizeType
Definition: element.h:94
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
std::size_t IndexType
Definition: flags.h:74
This is a 2D-2node truss element with 2 translational dofs per node.
Definition: geo_truss_element_base.hpp:37
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(GeoTrussElementBase)
void CalculateGeometricStiffnessMatrix(BoundedMatrix< double, TDim *TNumNodes, TDim *TNumNodes > &rGeometricStiffnessMatrix, const ProcessInfo &rCurrentProcessInfo)
This function assembles the geometric stiffness part of the total stiffness matrix.
void CalculateElasticStiffnessMatrix(MatrixType &rElasticStiffnessMatrix, const ProcessInfo &rCurrentProcessInfo)
This function assembles the elastic stiffness part of the total stiffness matrix.
void CreateTransformationMatrix(BoundedMatrix< double, TDim *TNumNodes, TDim *TNumNodes > &rRotationMatrix)
This function calculates the transformation matrix to globalize vectors and/or matrices.
virtual void WriteTransformationCoordinates(FullDofVectorType &rReferenceCoordinates)
This function calculates the current nodal postion for the transformation matrix.
GeoTrussElementBase()
Definition: geo_truss_element_base.hpp:60
This defines the geometrical object, base definition of the element and condition entities.
Definition: geometrical_object.h:58
Geometry base class.
Definition: geometry.h:71
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
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
#define DIM
Definition: logging_settings.hpp:38
std::size_t IndexType
Definition: binary_expression.cpp:25
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
TDataType Calculate(GeometryType &dummy, const Variable< TDataType > &rVariable)
Definition: add_geometries_to_python.cpp:103
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
Properties PropertiesType
Definition: regenerate_pfem_pressure_conditions_process.h:26
Geometry< Node > GeometryType
The definition of the geometry.
Definition: mortar_classes.h:37
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
def load(f)
Definition: ode_solve.py:307