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.
transient_Pw_element.hpp
Go to the documentation of this file.
1 // KRATOS___
2 // // ) )
3 // // ___ ___
4 // // ____ //___) ) // ) )
5 // // / / // // / /
6 // ((____/ / ((____ ((___/ / MECHANICS
7 //
8 // License: geo_mechanics_application/license.txt
9 //
10 // Main authors: Vahid Galavi
11 //
12 
13 #if !defined(KRATOS_GEO_PW_ELEMENT_H_INCLUDED)
14 #define KRATOS_GEO_PW_ELEMENT_H_INCLUDED
15 
16 // Project includes
17 #include "includes/serializer.h"
18 
19 // Application includes
20 #include "custom_elements/U_Pw_small_strain_element.hpp"
21 #include "custom_utilities/element_utilities.hpp"
24 
25 namespace Kratos
26 {
27 
28 template <unsigned int TDim, unsigned int TNumNodes>
29 class KRATOS_API(GEO_MECHANICS_APPLICATION) TransientPwElement : public UPwSmallStrainElement<TDim, TNumNodes>
30 {
31 public:
33 
35 
36  using IndexType = std::size_t;
38  using NodeType = Node;
41  using VectorType = Vector;
42  using MatrixType = Matrix;
45 
47  using SizeType = std::size_t;
48  using BaseType::CalculateRetentionResponse;
49  using BaseType::mConstitutiveLawVector;
50  using BaseType::mIsInitialised;
51  using BaseType::mRetentionLawVector;
52 
53  using ElementVariables = typename BaseType::ElementVariables;
54 
56 
58  TransientPwElement(IndexType NewId = 0) : BaseType(NewId) {}
59 
61  TransientPwElement(IndexType NewId, const NodesArrayType& ThisNodes)
62  : BaseType(NewId, ThisNodes)
63  {
64  }
65 
67  TransientPwElement(IndexType NewId, GeometryType::Pointer pGeometry)
68  : BaseType(NewId, pGeometry)
69  {
70  }
71 
73  TransientPwElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
74  : BaseType(NewId, pGeometry, pProperties)
75  {
76  }
77 
79  ~TransientPwElement() override {}
80 
82 
83  Element::Pointer Create(IndexType NewId,
84  NodesArrayType const& ThisNodes,
85  PropertiesType::Pointer pProperties) const override;
86 
87  Element::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override;
88 
90  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
91 
92  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
93 
94  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
95 
96  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
97 
98  void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
99 
100  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
101 
102  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
103 
104  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
105 
106  void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) override;
107 
108  void CalculateDampingMatrix(MatrixType& rDampingMatrix, const ProcessInfo& rCurrentProcessInfo) override;
109 
110  void GetValuesVector(Vector& rValues, int Step = 0) const override;
111 
112  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
113 
114  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
115 
117 
118  void CalculateOnIntegrationPoints(const Variable<double>& rVariable,
119  std::vector<double>& rOutput,
120  const ProcessInfo& rCurrentProcessInfo) override;
121 
123  std::vector<array_1d<double, 3>>& rOutput,
124  const ProcessInfo& rCurrentProcessInfo) override;
125 
126  void CalculateOnIntegrationPoints(const Variable<Matrix>& rVariable,
127  std::vector<Matrix>& rOutput,
128  const ProcessInfo& rCurrentProcessInfo) override;
129 
130  // Turn back information as a string.
131  std::string Info() const override
132  {
133  std::stringstream buffer;
134  buffer << "transient Pw flow Element #" << this->Id()
135  << "\nRetention law: " << mRetentionLawVector[0]->Info();
136  return buffer.str();
137  }
138 
139  // Print information about this object.
140  void PrintInfo(std::ostream& rOStream) const override
141  {
142  rOStream << "transient Pw flow Element #" << this->Id()
143  << "\nRetention law: " << mRetentionLawVector[0]->Info();
144  }
145 
147 
148 protected:
150 
152 
153  void CalculateAll(MatrixType& rLeftHandSideMatrix,
154  VectorType& rRightHandSideVector,
155  const ProcessInfo& CurrentProcessInfo,
156  const bool CalculateStiffnessMatrixFlag,
157  const bool CalculateResidualVectorFlag) override;
158 
159  void InitializeElementVariables(ElementVariables& rVariables, const ProcessInfo& CurrentProcessInfo) override;
160 
161  void CalculateAndAddLHS(MatrixType& rLeftHandSideMatrix, ElementVariables& rVariables) override;
162 
163  void CalculateAndAddRHS(VectorType& rRightHandSideVector, ElementVariables& rVariables, unsigned int GPoint) override;
164 
165  void CalculateKinematics(ElementVariables& rVariables, unsigned int PointNumber) override;
166 
167  void CalculateAndAddCompressibilityMatrix(MatrixType& rLeftHandSideMatrix, ElementVariables& rVariables) override;
168  void CalculateAndAddPermeabilityFlow(VectorType& rRightHandSideVector, ElementVariables& rVariables) override;
169  void CalculateAndAddFluidBodyFlow(VectorType& rRightHandSideVector, ElementVariables& rVariables) override;
170  void CalculateAndAddPermeabilityMatrix(MatrixType& rLeftHandSideMatrix, ElementVariables& rVariables) override;
171  void CalculateAndAddCompressibilityFlow(VectorType& rRightHandSideVector, ElementVariables& rVariables) override;
172 
173  unsigned int GetNumberOfDOF() const override;
175 
176 private:
178 
180 
182 
183  friend class Serializer;
184 
185  void save(Serializer& rSerializer) const override
186  {
188  }
189 
190  void load(Serializer& rSerializer) override{KRATOS_SERIALIZE_LOAD_BASE_CLASS(rSerializer, Element)}
191 
192  // Assignment operator.
193  TransientPwElement&
194  operator=(TransientPwElement const& rOther);
195 
196  // Copy constructor.
197  TransientPwElement(TransientPwElement const& rOther);
198 
199 }; // Class TransientPwElement
200 
201 } // namespace Kratos
202 
203 #endif // KRATOS_GEO_PW_ELEMENT_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Base class for all Elements.
Definition: element.h:60
std::size_t SizeType
Definition: element.h:94
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
This defines the geometrical object, base definition of the element and condition entities.
Definition: geometrical_object.h:58
Geometry base class.
Definition: geometry.h:71
PointerVector< TPointType > PointsArrayType
Definition: geometry.h:118
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
Definition: transient_Pw_element.hpp:30
TransientPwElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Constructor using Properties.
Definition: transient_Pw_element.hpp:73
TransientPwElement(IndexType NewId, const NodesArrayType &ThisNodes)
Constructor using an array of nodes.
Definition: transient_Pw_element.hpp:61
typename BaseType::ElementVariables ElementVariables
Definition: transient_Pw_element.hpp:53
~TransientPwElement() override
Destructor.
Definition: transient_Pw_element.hpp:79
TransientPwElement(IndexType NewId=0)
Default Constructor.
Definition: transient_Pw_element.hpp:58
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: transient_Pw_element.hpp:140
TransientPwElement(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor using Geometry.
Definition: transient_Pw_element.hpp:67
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(TransientPwElement)
std::string Info() const override
Turn back information as a string.
Definition: transient_Pw_element.hpp:131
Definition: U_Pw_small_strain_element.hpp:30
#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
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
pybind11::list CalculateOnIntegrationPoints(TObject &dummy, const Variable< TDataType > &rVariable, const ProcessInfo &rProcessInfo)
Definition: add_mesh_to_python.cpp:142
void InitializeSolutionStep(ConstructionUtility &rThisUtil, std::string ThermalSubModelPartName, std::string MechanicalSubModelPartName, std::string HeatFluxSubModelPartName, std::string HydraulicPressureSubModelPartName, bool thermal_conditions, bool mechanical_conditions, int phase)
Definition: add_custom_utilities_to_python.cpp:45
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
Internals::Matrix< double, AMatrix::dynamic, 1 > Vector
Definition: amatrix_interface.h:472
Internals::Matrix< double, AMatrix::dynamic, AMatrix::dynamic > Matrix
Definition: amatrix_interface.h:470
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
def load(f)
Definition: ode_solve.py:307