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.
UP_condition.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: Lorenzo Gracia
11 //
12 
13 
14 #if !defined(KRATOS_UP_CONDITION_H_INCLUDED )
15 #define KRATOS_UP_CONDITION_H_INCLUDED
16 
17 // System includes
18 #include <cmath>
19 
20 // Project includes
21 #include "includes/define.h"
22 #include "includes/condition.h"
23 #include "includes/serializer.h"
24 #include "includes/process_info.h"
25 
26 // Application includes
29 
30 namespace Kratos
31 {
32 
33 template< unsigned int TDim, unsigned int TNumNodes >
34 class UPCondition : public Condition
35 {
36 
37 public:
38 
40 
41 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
42 
43  // Default constructor
45 
46  // Constructor 1
47  UPCondition( IndexType NewId, GeometryType::Pointer pGeometry ) : Condition(NewId, pGeometry) {}
48 
49  // Constructor 2
50  UPCondition( IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties ) : Condition(NewId, pGeometry, pProperties)
51  {
53  }
54 
55  // Destructor
56  virtual ~UPCondition() {}
57 
58 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
59 
60  Condition::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
61 
62  void GetDofList(DofsVectorType& rConditionDofList, const ProcessInfo& rCurrentProcessInfo) const override;
63 
64 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
65 
66  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
67 
68  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo) override;
69 
70  void CalculateRightHandSide(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
71 
72  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
73 
74 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
75 
76 protected:
77 
79  {
81  double Density;
82 
85 
91 
95 
101  };
102 
103  // Member Variables
104 
106 
107 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
108 
109  void CalculateAll( MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo);
110 
111  void CalculateLHS( MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo);
112 
113  void CalculateRHS( VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo);
114 
115  void CalculateNormalVector(VectorType& rNormalVector, const Matrix& Jacobian);
116 
117  void CalculateIntegrationCoefficient(double& rIntegrationCoefficient, const Matrix& Jacobian, const double& weight);
118 
119  void InitializeConditionVariables(ConditionVariables& rVariables, const GeometryType& Geom, const PropertiesType& Prop, const ProcessInfo& rCurrentProcessInfo);
120 
121  void GetAccelerationVector( Vector& rValues, int Step );
122 
123  void CalculateLHSContribution(MatrixType& rLeftHandSideMatrix, ConditionVariables& rVariables);
124 
125  void CalculateRHSContribution(VectorType& rRightHandSideVector, ConditionVariables& rVariables);
126 
127 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
128 
129 private:
130 
131  // Serialization
132 
133  friend class Serializer;
134 
135  void save(Serializer& rSerializer) const override
136  {
138  }
139 
140  void load(Serializer& rSerializer) override
141  {
143  }
144 
145 }; // class UPCondition.
146 
147 } // namespace Kratos.
148 
149 #endif // KRATOS_UP_CONDITION_H_INCLUDED defined
Base class for all Conditions.
Definition: condition.h:59
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
std::vector< DofType::Pointer > DofsVectorType
Definition: condition.h:100
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
Geometry base class.
Definition: geometry.h:71
IntegrationMethod GetDefaultIntegrationMethod() const
Definition: geometry.h:2004
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
Definition: UP_condition.hpp:35
void CalculateIntegrationCoefficient(double &rIntegrationCoefficient, const Matrix &Jacobian, const double &weight)
void CalculateLHSContribution(MatrixType &rLeftHandSideMatrix, ConditionVariables &rVariables)
Definition: UP_condition.cpp:501
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
void CalculateAll(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: UP_condition.cpp:233
void GetAccelerationVector(Vector &rValues, int Step)
Definition: UP_condition.cpp:482
KRATOS_CLASS_POINTER_DEFINITION(UPCondition)
virtual ~UPCondition()
Definition: UP_condition.hpp:56
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
It creates a new condition pointer.
Definition: UP_condition.cpp:21
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: UP_condition.cpp:143
void CalculateRHS(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: UP_condition.cpp:293
void GetDofList(DofsVectorType &rConditionDofList, const ProcessInfo &rCurrentProcessInfo) const override
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: UP_condition.cpp:100
void CalculateRHSContribution(VectorType &rRightHandSideVector, ConditionVariables &rVariables)
Definition: UP_condition.cpp:520
void CalculateLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: UP_condition.cpp:247
UPCondition()
Definition: UP_condition.hpp:44
void InitializeConditionVariables(ConditionVariables &rVariables, const GeometryType &Geom, const PropertiesType &Prop, const ProcessInfo &rCurrentProcessInfo)
Definition: UP_condition.cpp:451
void CalculateNormalVector(VectorType &rNormalVector, const Matrix &Jacobian)
UPCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Definition: UP_condition.hpp:50
void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: UP_condition.cpp:124
GeometryData::IntegrationMethod mThisIntegrationMethod
Definition: UP_condition.hpp:105
UPCondition(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: UP_condition.hpp:47
#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
Definition: UP_condition.hpp:79
array_1d< double, TNumNodes *TDim > UVector
Definition: UP_condition.hpp:99
Vector Np
Variables computed at each GP.
Definition: UP_condition.hpp:87
double AcelerationCoefficient
ProcessInfo variables.
Definition: UP_condition.hpp:84
array_1d< double, TNumNodes > PVector
Definition: UP_condition.hpp:100
BoundedMatrix< double, TDim, TNumNodes *TDim > Nu
Definition: UP_condition.hpp:89
double IntegrationCoefficient
Definition: UP_condition.hpp:90
Vector AccelerationVector
Definition: UP_condition.hpp:94
BoundedMatrix< double, TNumNodes *TDim, TNumNodes > UPMatrix
Auxiliary Variables.
Definition: UP_condition.hpp:97
Vector NormalVector
Definition: UP_condition.hpp:88
array_1d< double, TNumNodes > PressureVector
Nodal variables.
Definition: UP_condition.hpp:93
BoundedMatrix< double, TNumNodes, TNumNodes *TDim > PUMatrix
Definition: UP_condition.hpp:98
double Density
Properties variables.
Definition: UP_condition.hpp:81