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.
mass_element.h
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: Philipp Bucher (https://github.com/philbucher)
10 //
11 
12 #pragma once
13 
14 // System includes
15 
16 // External includes
17 
18 // Project includes
19 #include "includes/element.h"
20 
21 namespace Kratos
22 {
25 
26 class MassElement : public Element
27 {
28 public:
29 
32 
34 
37 
40 
44 
49  : Element(NewId) {}
50 
54  MassElement(IndexType NewId, const NodesArrayType& ThisNodes)
55  : Element(NewId, ThisNodes) {}
56 
60  MassElement(IndexType NewId, GeometryType::Pointer pGeometry)
61  : Element(NewId, pGeometry) {}
62 
66  MassElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
67  : Element(NewId, pGeometry, pProperties) {}
68 
72  MassElement(MassElement const& rOther)
73  : Element(rOther) {}
74 
78  ~MassElement() override = default;
79 
83 
84 
88 
96  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
97 
105  Element::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override;
106 
114  Element::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override;
115 
121  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& CurrentProcessInfo) const override;
122 
128  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& CurrentProcessInfo) const override;
129 
130 
134  void GetValuesVector(Vector& values, int Step = 0) const override;
135 
139  void GetFirstDerivativesVector(Vector& values, int Step = 0) const override;
140 
144  void GetSecondDerivativesVector(Vector& values, int Step = 0) const override;
145 
151  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
152 
162  MatrixType& rLeftHandSideMatrix,
163  VectorType& rRightHandSideVector,
164  const ProcessInfo& rCurrentProcessInfo) override;
165 
172  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo) override;
173 
180  void CalculateRightHandSide(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
181 
188  void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) override;
189 
196  void CalculateDampingMatrix(MatrixType& rDampingMatrix, const ProcessInfo& rCurrentProcessInfo) override;
197 
207  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
208 
212 
214  std::string Info() const override;
215 
217  void PrintInfo(std::ostream& rOStream) const override;
218 
220  void PrintData(std::ostream& rOStream) const override;
221 
223 
224 private:
225 
228 
229  double mMass = 0.0;
230 
232 
235 
236  void GenericGetValuesVector(Vector& rValues, int Step, const ArrayVariableType& rVariable) const;
237 
238  double GetElementMass() const;
239 
243 
244  friend class Serializer;
245 
246  void save(Serializer& rSerializer) const override;
247  void load(Serializer& rSerializer) override;
248 
250 
251 }; // Class MassElement
252 
254 
255 } // namespace Kratos.
Base class for all Elements.
Definition: element.h:60
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
Definition: mass_element.h:27
MassElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Definition: mass_element.h:66
void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &CurrentProcessInfo) const override
Definition: mass_element.cpp:81
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: mass_element.cpp:149
Element::Pointer Clone(IndexType NewId, NodesArrayType const &ThisNodes) const override
Definition: mass_element.cpp:52
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(MassElement)
Pointer definition of MassElement.
void GetSecondDerivativesVector(Vector &values, int Step=0) const override
Definition: mass_element.cpp:135
MassElement(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: mass_element.h:60
Variable< array_1d< double, 3 > > ArrayVariableType
Definition: mass_element.h:33
MassElement(MassElement const &rOther)
Definition: mass_element.h:72
void GetValuesVector(Vector &values, int Step=0) const override
Definition: mass_element.cpp:125
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
Definition: mass_element.cpp:140
MassElement(IndexType NewId, const NodesArrayType &ThisNodes)
Definition: mass_element.h:54
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &CurrentProcessInfo) const override
Definition: mass_element.cpp:59
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: mass_element.cpp:321
void GetFirstDerivativesVector(Vector &values, int Step=0) const override
Definition: mass_element.cpp:130
MassElement(IndexType NewId=0)
Definition: mass_element.h:48
void CalculateDampingMatrix(MatrixType &rDampingMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: mass_element.cpp:225
~MassElement() override=default
void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: mass_element.cpp:195
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: mass_element.cpp:168
std::string Info() const override
Turn back information as a string.
Definition: mass_element.cpp:312
void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: mass_element.cpp:158
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Definition: mass_element.cpp:32
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Definition: mass_element.cpp:262
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: mass_element.cpp:328
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
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
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
list values
Definition: bombardelli_test.py:42
def load(f)
Definition: ode_solve.py:307