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.
free_surface_condition.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosDamApplication $
3 // Last Modified by: $Author: Lorenzo Gracia $
4 // Date: $Date: January 2016 $
5 // Revision: $Revision: 1.0 $
6 //
7 
8 #if !defined(KRATOS_FREE_SURFACE_CONDITION_H_INCLUDED )
9 #define KRATOS_FREE_SURFACE_CONDITION_H_INCLUDED
10 
11 // System includes
12 #include <cmath>
13 
14 // Project includes
15 #include "includes/define.h"
16 #include "includes/condition.h"
17 #include "includes/serializer.h"
18 #include "includes/process_info.h"
19 
20 // Application includes
22 
23 namespace Kratos
24 {
25 
26 template< unsigned int TDim, unsigned int TNumNodes >
28 {
29 
30 public:
31 
33 
34  typedef std::size_t IndexType;
36  typedef Node NodeType;
39  typedef Vector VectorType;
40  typedef Matrix MatrixType;
41 
42 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
43 
44  // Default constructor
46 
47  // Constructor 1
48  FreeSurfaceCondition( IndexType NewId, GeometryType::Pointer pGeometry ) : Condition(NewId, pGeometry) {}
49 
50  // Constructor 2
51  FreeSurfaceCondition( IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties ) : Condition(NewId, pGeometry, pProperties)
52  {
54  }
55 
56  // Destructor
57  virtual ~FreeSurfaceCondition() {}
58 
59 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
60 
61  Condition::Pointer Create(IndexType NewId,NodesArrayType const& ThisNodes,PropertiesType::Pointer pProperties ) const override;
62 
63  void GetDofList(DofsVectorType& rConditionDofList, const ProcessInfo& rCurrentProcessInfo ) const override;
64 
65 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
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 EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo ) const override;
74 
75 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
76 
77 protected:
78 
79  // Member Variables
80 
82 
83 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
84 
85  void CalculateAll( MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo );
86 
87  virtual void CalculateLHS( MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo );
88 
89  virtual void CalculateRHS( VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo );
90 
91  void CalculateIntegrationCoefficient(double& rIntegrationCoefficient, const Matrix& Jacobian, const double& weight);
92 
93 //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
94 
95 private:
96 
97  // Serialization
98 
99  friend class Serializer;
100 
101  void save(Serializer& rSerializer) const override
102  {
104  }
105 
106  void load(Serializer& rSerializer) override
107  {
109  }
110 
111 }; // class FreeSurfaceCondition.
112 
113 } // namespace Kratos.
114 
115 #endif // KRATOS_FREE_SURFACE_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
Definition: free_surface_condition.hpp:28
GeometryData::IntegrationMethod mThisIntegrationMethod
Definition: free_surface_condition.hpp:81
FreeSurfaceCondition(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: free_surface_condition.hpp:48
Vector VectorType
Definition: free_surface_condition.hpp:39
void CalculateAll(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: free_surface_condition.cpp:127
Node NodeType
Definition: free_surface_condition.hpp:36
void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Definition: free_surface_condition.cpp:68
KRATOS_CLASS_POINTER_DEFINITION(FreeSurfaceCondition)
FreeSurfaceCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Definition: free_surface_condition.hpp:51
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Definition: free_surface_condition.cpp:106
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition: free_surface_condition.hpp:38
virtual void CalculateLHS(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: free_surface_condition.cpp:137
FreeSurfaceCondition()
Definition: free_surface_condition.hpp:45
std::size_t IndexType
Definition: free_surface_condition.hpp:34
virtual void CalculateRHS(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: free_surface_condition.cpp:182
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
It creates a new condition pointer.
Definition: free_surface_condition.cpp:15
Matrix MatrixType
Definition: free_surface_condition.hpp:40
void CalculateIntegrationCoefficient(double &rIntegrationCoefficient, const Matrix &Jacobian, const double &weight)
void GetDofList(DofsVectorType &rConditionDofList, const ProcessInfo &rCurrentProcessInfo) const override
Definition: free_surface_condition.cpp:23
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: free_surface_condition.cpp:44
Geometry< NodeType > GeometryType
Definition: free_surface_condition.hpp:37
virtual ~FreeSurfaceCondition()
Definition: free_surface_condition.hpp:57
Properties PropertiesType
Definition: free_surface_condition.hpp:35
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Definition: free_surface_condition.cpp:87
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
This class defines the node.
Definition: node.h:65
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
#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