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.
U_Pw_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: Ignasi de Pouplana
11 //
12 
13 
14 #if !defined(KRATOS_U_PW_CONDITION_H_INCLUDED )
15 #define KRATOS_U_PW_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 KRATOS_API(POROMECHANICS_APPLICATION) UPwCondition : public Condition
35 {
36 
37 public:
38 
41 
44 
46 
47 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
48 
49  // Default constructor
51 
52  // Constructor 1
53  UPwCondition( IndexType NewId, GeometryType::Pointer pGeometry ) : Condition(NewId, pGeometry) {}
54 
55  // Constructor 2
56  UPwCondition( IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties ) : Condition(NewId, pGeometry, pProperties)
57  {
58  mThisIntegrationMethod = this->GetIntegrationMethod();
59  }
60 
61  // Destructor
62  virtual ~UPwCondition() {}
63 
64 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
65 
66  Condition::Pointer Create(IndexType NewId,NodesArrayType const& ThisNodes,PropertiesType::Pointer pProperties ) const override;
67 
68  void GetDofList(DofsVectorType& rConditionDofList,const ProcessInfo& rCurrentProcessInfo ) const override;
69 
70 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
71 
72  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,VectorType& rRightHandSideVector,const ProcessInfo& rCurrentProcessInfo ) override;
73 
74  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix,const ProcessInfo& rCurrentProcessInfo ) override;
75 
76  void CalculateRightHandSide(VectorType& rRightHandSideVector,const ProcessInfo& rCurrentProcessInfo ) override;
77 
78  void EquationIdVector(EquationIdVectorType& rResult,const ProcessInfo& rCurrentProcessInfo ) const override;
79 
80 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
81 
82  void AddExplicitContribution(const VectorType& rRHSVector,
83  const Variable<VectorType>& rRHSVariable,
84  const Variable<array_1d<double,3> >& rDestinationVariable,
85  const ProcessInfo& rCurrentProcessInfo
86  ) override;
87 
88 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
89 
90 protected:
91 
92  // Member Variables
93 
94  GeometryData::IntegrationMethod mThisIntegrationMethod;
95 
96 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
97 
98  virtual void CalculateAll( MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo );
99 
100  virtual void CalculateRHS( VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo );
101 
102 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
103 
104 private:
105 
106  // Serialization
107 
108  friend class Serializer;
109 
110  void save(Serializer& rSerializer) const override
111  {
113  }
114 
115  void load(Serializer& rSerializer) override
116  {
118  }
119 
120 }; // class UPwCondition.
121 
122 } // namespace Kratos.
123 
124 #endif // KRATOS_U_PW_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
IntegrationMethod
Definition: geometry_data.h:76
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: U_Pw_condition.hpp:37
void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
virtual void CalculateRHS(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
UPwCondition(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: U_Pw_condition.hpp:53
UPwCondition()
Definition: U_Pw_condition.hpp:50
UPwCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Definition: U_Pw_condition.hpp:56
BaseType::SizeType SizeType
Definition of the size type.
Definition: U_Pw_condition.hpp:43
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
It creates a new condition pointer.
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Condition BaseType
We define the base class Condition.
Definition: U_Pw_condition.hpp:40
KRATOS_CLASS_POINTER_DEFINITION(UPwCondition)
virtual void CalculateAll(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
void GetDofList(DofsVectorType &rConditionDofList, const ProcessInfo &rCurrentProcessInfo) const override
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
virtual ~UPwCondition()
Definition: U_Pw_condition.hpp:62
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#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 SizeType
Definition: nurbs_utilities.h:41
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
def load(f)
Definition: ode_solve.py:307