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.
wave_equation_element.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosDamApplication $
3 // Last modified by: $Author: lgracia $
4 // Date: $Date: December 2016 $
5 // Revision: $Revision: 1.0 $
6 //
7 //
8 
9 #if !defined(KRATOS_WAVE_EQUATION_ELEMENT_H_INCLUDED )
10 #define KRATOS_WAVE_EQUATION_ELEMENT_H_INCLUDED
11 
12 // Project includes
13 #include "containers/array_1d.h"
14 #include "includes/define.h"
15 #include "includes/element.h"
16 #include "includes/serializer.h"
17 #include "geometries/geometry.h"
18 #include "utilities/math_utils.h"
19 
21 
22 namespace Kratos
23 {
24 
25 template< unsigned int TDim, unsigned int TNumNodes >
27 {
28 
29 public:
30 
32 
34 
36  WaveEquationElement(IndexType NewId = 0) : Element( NewId ) {}
37 
39  WaveEquationElement(IndexType NewId, const NodesArrayType& ThisNodes) : Element(NewId, ThisNodes) {}
40 
42  WaveEquationElement(IndexType NewId, GeometryType::Pointer pGeometry) : Element( NewId, pGeometry ) {}
43 
45  WaveEquationElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties) : Element( NewId, pGeometry, pProperties )
46  {
48  }
49 
51  virtual ~WaveEquationElement() {}
52 
54 
55  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
56 
57  Element::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override;
58 
59  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
60 
61  void GetDofList( DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
62 
63  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
64 
66 
67  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo ) override;
68 
69  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo ) override;
70 
71  void CalculateRightHandSide(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo ) override;
72 
73  void GetValuesVector(Vector& rValues, int Step = 0) const override;
74 
75  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
76 
77  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
78 
80 
81 protected:
82 
84 
86 
88 
89  void CalculateAll(MatrixType& rLeftHandSideMatrix,VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo );
90 
91  void CalculateLHS(MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo );
92 
93  void CalculateRHS(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo );
94 
95  void CalculateIntegrationCoefficient(double& rIntegrationCoefficient, const double& detJ, const double& weight);
96 
98 
99 private:
100 
102 
103  friend class Serializer;
104 
105  void save(Serializer& rSerializer) const override
106  {
108  }
109 
110  void load(Serializer& rSerializer) override
111  {
113  }
114 
117 
120 
121 
122 }; // Class WaveEquationElement
123 
124 } // namespace Kratos
125 
126 #endif // KRATOS_WAVE_EQUATION_ELEMENT_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Base class for all Elements.
Definition: element.h:60
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
Matrix MatrixType
Definition: element.h:90
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
std::size_t IndexType
Definition: flags.h:74
GeometryType & GetGeometry()
Returns the reference of the geometry.
Definition: geometrical_object.h:158
IntegrationMethod
Definition: geometry_data.h:76
IntegrationMethod GetDefaultIntegrationMethod() const
Definition: geometry.h:2004
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
Definition: wave_equation_element.hpp:27
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Definition: wave_equation_element.cpp:98
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Definition: wave_equation_element.cpp:32
WaveEquationElement(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor using Geometry.
Definition: wave_equation_element.hpp:42
void GetFirstDerivativesVector(Vector &rValues, int Step=0) const override
Definition: wave_equation_element.cpp:201
GeometryData::IntegrationMethod mThisIntegrationMethod
Member Variables.
Definition: wave_equation_element.hpp:85
virtual ~WaveEquationElement()
Destructor.
Definition: wave_equation_element.hpp:51
void CalculateAll(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: wave_equation_element.cpp:238
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Definition: wave_equation_element.cpp:16
WaveEquationElement(IndexType NewId, const NodesArrayType &ThisNodes)
Constructor using an array of nodes.
Definition: wave_equation_element.hpp:39
void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: wave_equation_element.cpp:144
void GetValuesVector(Vector &rValues, int Step=0) const override
Definition: wave_equation_element.cpp:182
WaveEquationElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Constructor using Properties.
Definition: wave_equation_element.hpp:45
void CalculateRHS(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: wave_equation_element.cpp:297
void GetSecondDerivativesVector(Vector &rValues, int Step=0) const override
Definition: wave_equation_element.cpp:220
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: wave_equation_element.cpp:162
void CalculateIntegrationCoefficient(double &rIntegrationCoefficient, const double &detJ, const double &weight)
Definition: wave_equation_element.cpp:355
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(WaveEquationElement)
void CalculateLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: wave_equation_element.cpp:250
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: wave_equation_element.cpp:119
void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &rCurrentProcessInfo) const override
Definition: wave_equation_element.cpp:77
WaveEquationElement(IndexType NewId=0)
Default Constructor.
Definition: wave_equation_element.hpp:36
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
def load(f)
Definition: ode_solve.py:307