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.
updated_lagrangian_element.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosSolidMechanicsApplication $
3 // Created by: $Author: JMCarbonell $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: July 2013 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_UPDATED_LAGRANGIAN_ELEMENT_H_INCLUDED)
11 #define KRATOS_UPDATED_LAGRANGIAN_ELEMENT_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
19 
20 
21 namespace Kratos
22 {
37 
39 
45 class KRATOS_API(SOLID_MECHANICS_APPLICATION) UpdatedLagrangianElement
46  : public LargeDisplacementElement
47 {
48 public:
49 
55  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
62 
68 
70  UpdatedLagrangianElement(IndexType NewId, GeometryType::Pointer pGeometry);
71 
72  UpdatedLagrangianElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
73 
76 
79 
83 
86 
90 
101  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
102 
103 
111  Element::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override;
112 
113 
114  // //************* GETTING METHODS
115 
116  //SET
117 
121  void SetValuesOnIntegrationPoints(const Variable<double>& rVariable, const std::vector<double>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
122 
123 
124  //GET:
125 
129  void CalculateOnIntegrationPoints(const Variable<double>& rVariable, std::vector<double>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
130 
131 
132 
133  //************* STARTING - ENDING METHODS
134 
139  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
140 
141  //************************************************************************************
142  //************************************************************************************
150  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
151 
155 
163  std::string Info() const override
164  {
165  std::stringstream buffer;
166  buffer << "Updated Lagrangian Element #" << Id();
167  return buffer.str();
168  }
169 
171  void PrintInfo(std::ostream& rOStream) const override
172  {
173  rOStream << "Updated Lagrangian Element #" << Id();
174  }
175 
177  void PrintData(std::ostream& rOStream) const override
178  {
179  GetGeometry().PrintData(rOStream);
180  }
181 
186 
187 protected:
193 
197  std::vector< Matrix > mDeformationGradientF0;
198 
203 
208  {
209  }
210 
214  void InitializeElementData(ElementDataType& rVariables,
215  const ProcessInfo& rCurrentProcessInfo) override;
216 
217 
218 
222  void FinalizeStepVariables(ElementDataType & rVariables,
223  const double& rPointNumber ) override;
224 
225 
229  void CalculateKinematics(ElementDataType& rVariables,
230  const double& rPointNumber) override;
231 
232 
236  void CalculateKinetics(ElementDataType& rVariables,
237  const double& rPointNumber) override;
238 
239 
240 
244  void CalculateDeformationGradient(Matrix& rF,
245  const Matrix& rDN_DX,
246  const Matrix& rDeltaPosition);
247 
248 
252  void GetHistoricalVariables( ElementDataType& rVariables,
253  const double& rPointNumber ) override;
254 
258  double& CalculateVolumeChange(double& rVolumeChange, ElementDataType& rVariables) override;
259 
270 
271 private:
272 
288 
292  friend class Serializer;
293 
294  void save(Serializer& rSerializer) const override;
295 
296  void load(Serializer& rSerializer) override;
297 
304 
305 }; // Class UpdatedLagrangianElement
306 
314 
315 } // namespace Kratos.
316 #endif // KRATOS_UPDATED_LAGRANGIAN_ELEMENT_H_INCLUDED defined
Definition: constitutive_law.h:47
std::size_t IndexType
Definition: flags.h:74
IntegrationMethod
Definition: geometry_data.h:76
std::size_t SizeType
Definition: geometry_data.h:173
Large Displacement Lagrangian Element for 3D and 2D geometries. (base class)
Definition: large_displacement_element.hpp:46
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
A stabilized element for the incompressible Navier-Stokes equations.
Definition: updated_lagrangian_element.h:64
UpdatedLagrangianElement()
Definition: updated_lagrangian_element.hpp:207
UpdatedLagrangianElement & operator=(UpdatedLagrangianElement const &rOther)
Assignment operator.
void Initialize(const ProcessInfo &rCurrentProcessInfo) override
std::string Info() const override
Turn back information as a string.
Definition: updated_lagrangian_element.hpp:163
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(UpdatedLagrangianElement)
Counted pointer of UpdatedLagrangianElement.
~UpdatedLagrangianElement() override
Destructor.
ConstitutiveLaw ConstitutiveLawType
Definition: updated_lagrangian_element.hpp:53
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: updated_lagrangian_element.hpp:171
LargeDisplacementElement::ElementDataType ElementDataType
Type for element variables.
Definition: updated_lagrangian_element.hpp:61
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: updated_lagrangian_element.hpp:57
UpdatedLagrangianElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
UpdatedLagrangianElement(UpdatedLagrangianElement const &rOther)
Copy constructor.
std::vector< Matrix > mDeformationGradientF0
Definition: updated_lagrangian_element.hpp:197
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: updated_lagrangian_element.hpp:55
GeometryData::SizeType SizeType
Type for size.
Definition: updated_lagrangian_element.hpp:59
UpdatedLagrangianElement(IndexType NewId, GeometryType::Pointer pGeometry)
Default constructors.
Vector mDeterminantF0
Definition: updated_lagrangian_element.hpp:202
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: updated_lagrangian_element.hpp:177
Element::Pointer Clone(IndexType NewId, NodesArrayType const &ThisNodes) const override
int Check(const ProcessInfo &rCurrentProcessInfo) const override
void SetValuesOnIntegrationPoints(TObject &dummy, const Variable< TDataType > &rVariable, const std::vector< TDataType > &values, const ProcessInfo &rCurrentProcessInfo)
Definition: add_mesh_to_python.cpp:185
pybind11::list CalculateOnIntegrationPoints(TObject &dummy, const Variable< TDataType > &rVariable, const ProcessInfo &rProcessInfo)
Definition: add_mesh_to_python.cpp:142
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
def load(f)
Definition: ode_solve.py:307
Definition: solid_element.hpp:83