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.
axisym_updated_lagrangian_U_wP_element.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosPfemSolidMechanicsApplication $
3 // Created by: $Author: LMonforte $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: July 2015 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_AXISYM_UPDATED_LAGRANGIAN_U_wP_ELEMENT_H_INCLUDED )
11 #define KRATOS_AXISYM_UPDATED_LAGRANGIAN_U_wP_ELEMENT_H_INCLUDED
12 
13 
14 // System includes
15 
16 // External includes
17 
18 // Project includes
20 
21 namespace Kratos
22 {
37 
39 
40 
41 class KRATOS_API(PFEM_SOLID_MECHANICS_APPLICATION) AxisymUpdatedLagrangianUwPElement
43 {
44 public:
45 
51  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
60 
64 
67 
70 
72  AxisymUpdatedLagrangianUwPElement(IndexType NewId, GeometryType::Pointer pGeometry);
73 
74  AxisymUpdatedLagrangianUwPElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
75 
78 
79 
82 
86 
89 
90 
94 
105  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
106 
114  Element::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override;
115 
116  //************* GETTING METHODS
117  //************************************************************************************
118  //************************************************************************************
126  int Check(const ProcessInfo& rCurrentProcessInfo) override;
127 
128  //GET:
129 
133  void GetValueOnIntegrationPoints(const Variable<double>& rVariable, std::vector<double>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
134 
135  void GetValueOnIntegrationPoints(const Variable<Vector>& rVariable, std::vector<Vector>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
136 
137  void GetValueOnIntegrationPoints( const Variable<Matrix>& rVariable, std::vector<Matrix>& rValue, const ProcessInfo& rCurrentProcessInfo) override;
138 
139  //************* STARTING - ENDING METHODS
140 
141 
145  void GetDofList(DofsVectorType& rElementalDofList, ProcessInfo& rCurrentProcessInfo) override;
146 
150  void EquationIdVector(EquationIdVectorType& rResult, ProcessInfo& rCurrentProcessInfo) override;
151 
155  void GetValuesVector(Vector& rValues, int Step = 0) const override;
156 
160  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
161 
165  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
166 
167 
168 
172 
183 protected:
189 
190 
191  /****
192  the time step (requiered). It shall be somewhere else.
193  ****/
194  double mTimeStep;
195 
199 
203 
204 
209  virtual void CalculateAndAddLHS(LocalSystemComponents& rLocalSystem,
210  ElementDataType& rVariables,
211  double& rIntegrationWeight) override;
212 
217  virtual void CalculateAndAddRHS(LocalSystemComponents& rLocalSystem,
218  ElementDataType& rVariables,
219  Vector& rVolumeForce,
220  double& rIntegrationWeight) override;
221 
225  virtual void InitializeElementData(ElementDataType & rVariables, const ProcessInfo& rCurrentProcessInfo) override;
226 
227 
228 
232  void InitializeSystemMatrices(MatrixType& rLeftHandSideMatrix,
233  VectorType& rRightHandSideVector,
234  Flags& rCalculationFlags) override;
235 
236  //on integration points:
240  void CalculateOnIntegrationPoints(const Variable<double>& rVariable, std::vector<double>& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
241 
242  void CalculateOnIntegrationPoints(const Variable<Vector>& rVariable, std::vector<Vector>& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
243 
244  void CalculateOnIntegrationPoints(const Variable<Matrix>& rVariable, std::vector<Matrix>& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
245 
256 
257 private:
258 
264 
265 
269 
270 
274 
275 
280 
284  friend class Serializer;
285 
286  // A private default constructor necessary for serialization
287 
288  virtual void save(Serializer& rSerializer) const override;
289 
290  virtual void load(Serializer& rSerializer) override;
291 
292 
299 
300 
301 }; // Class UpdatedLagrangianUwPElement
302 
303 
304 
305 } // namespace Kratos
306 #endif // KRATOS_UPDATED_LAGRANGIAN_U_wP_ELEMENT_H_INCLUDED
307 
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Axisymmetric Updated Lagrangian Large Displacement Lagrangian U-Pw Element.
Definition: axisym_updated_lagrangian_U_wP_element.hpp:43
ConstitutiveLawType::StressMeasure StressMeasureType
StressMeasure from constitutive laws.
Definition: axisym_updated_lagrangian_U_wP_element.hpp:53
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: axisym_updated_lagrangian_U_wP_element.hpp:51
ConstitutiveLaw ConstitutiveLawType
Definition: axisym_updated_lagrangian_U_wP_element.hpp:49
GeometryData::SizeType SizeType
Type for size.
Definition: axisym_updated_lagrangian_U_wP_element.hpp:57
AxisymmetricUpdatedLagrangianElement::ElementDataType ElementDataType
Type for element variables.
Definition: axisym_updated_lagrangian_U_wP_element.hpp:59
double mTimeStep
Definition: axisym_updated_lagrangian_U_wP_element.hpp:194
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: axisym_updated_lagrangian_U_wP_element.hpp:55
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(AxisymUpdatedLagrangianUwPElement)
Counted pointer of LargeDisplacementUPElement.
Axisymmetric Updated Lagrangian Element 2D geometries.
Definition: axisymmetric_updated_lagrangian_element.hpp:47
Definition: constitutive_law.h:47
StressMeasure
Definition: constitutive_law.h:69
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
Definition: flags.h:58
std::size_t IndexType
Definition: flags.h:74
IntegrationMethod
Definition: geometry_data.h:76
std::size_t SizeType
Definition: geometry_data.h:173
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
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
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
Definition: solid_element.hpp:233