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_U_P_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_U_P_ELEMENT_H_INCLUDED)
11 #define KRATOS_UPDATED_LAGRANGIAN_U_P_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) UpdatedLagrangianUPElement
47 {
48 public:
49 
55  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
64 
67 
71 
73  UpdatedLagrangianUPElement(IndexType NewId, GeometryType::Pointer pGeometry);
74 
75  UpdatedLagrangianUPElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
76 
79 
81  ~UpdatedLagrangianUPElement() override;
82 
86 
89 
93 
104  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
105 
113  Element::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override;
114 
115 
116  //************* GETTING METHODS
117 
118  //SET
119 
123  void SetValuesOnIntegrationPoints(const Variable<double>& rVariable, const std::vector<double>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
124 
125 
126  //GET:
127 
131  void CalculateOnIntegrationPoints(const Variable<double>& rVariable, std::vector<double>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
132 
133 
134 
135  //************* STARTING - ENDING METHODS
136 
141  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
142 
143  //************************************************************************************
144  //************************************************************************************
152  //int Check(const ProcessInfo& rCurrentProcessInfo);
153 
154 
158 
169 
170 protected:
176 
180  std::vector< Matrix > mDeformationGradientF0;
181 
186 
191  {
192  }
193 
197 
198 
203  void CalculateAndAddLHS(LocalSystemComponents& rLocalSystem,
204  ElementDataType& rVariables,
205  double& rIntegrationWeight) override;
206 
211  void CalculateAndAddRHS(LocalSystemComponents& rLocalSystem,
212  ElementDataType& rVariables,
213  Vector& rVolumeForce,
214  double& rIntegrationWeight) override;
215 
219  void InitializeElementData(ElementDataType & rVariables, const ProcessInfo& rCurrentProcessInfo) override;
220 
221 
225  void FinalizeStepVariables(ElementDataType & rVariables, const double& rPointNumber) override;
226 
227 
231  void CalculateKinematics(ElementDataType& rVariables,
232  const double& rPointNumber) override;
233 
234 
238  void CalculateDeformationGradient(Matrix& rF,
239  const Matrix& rDN_DX,
240  const Matrix& rDeltaPosition);
241 
242 
246  void GetHistoricalVariables( ElementDataType& rVariables,
247  const double& rPointNumber ) override;
248 
252  double& CalculateVolumeChange(double& rVolumeChange, ElementDataType& rVariables) override;
253 
264 
265 private:
266 
272 
273 
277 
278 
282 
283 
288 
292  friend class Serializer;
293 
294  // A private default constructor necessary for serialization
295 
296  void save(Serializer& rSerializer) const override;
297 
298  void load(Serializer& rSerializer) override;
299 
300 
307 
308 }; // Class UpdatedLagrangianUPElement
309 
317 
318 } // namespace Kratos.
319 #endif // KRATOS_UPDATED_LAGRANGIAN_U_P_ELEMENT_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Definition: constitutive_law.h:47
StressMeasure
Definition: constitutive_law.h:69
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 U-P Element for 3D and 2D geometries. Linear Triangles and Tetrahedra (...
Definition: large_displacement_U_P_element.hpp:47
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
Spatial Lagrangian U-P Element for 3D and 2D geometries. Linear Triangles and Tetrahedra.
Definition: updated_lagrangian_U_P_element.hpp:47
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: updated_lagrangian_U_P_element.hpp:59
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(UpdatedLagrangianUPElement)
Counted pointer of UpdatedLagrangianUPElement.
UpdatedLagrangianUPElement()
Definition: updated_lagrangian_U_P_element.hpp:190
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: updated_lagrangian_U_P_element.hpp:55
Vector mDeterminantF0
Definition: updated_lagrangian_U_P_element.hpp:185
LargeDisplacementElement::ElementDataType ElementDataType
Type for element variables.
Definition: updated_lagrangian_U_P_element.hpp:63
GeometryData::SizeType SizeType
Type for size.
Definition: updated_lagrangian_U_P_element.hpp:61
ConstitutiveLaw ConstitutiveLawType
Definition: updated_lagrangian_U_P_element.hpp:53
ConstitutiveLawType::StressMeasure StressMeasureType
StressMeasure from constitutive laws.
Definition: updated_lagrangian_U_P_element.hpp:57
std::vector< Matrix > mDeformationGradientF0
Definition: updated_lagrangian_U_P_element.hpp:180
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
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