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_J_element.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosSolidMechanicsApplication $
3 // Last modified by: $Author: LMonforte $
4 // Date: $Date: July 2015 $
5 // Revision: $Revision: 0.0 $
6 //
7 //
8 
9 #if !defined(KRATOS_AXISYM_UPDATED_LAGRANGIAN_U_J_ELEMENT_H_INCLUDED )
10 #define KRATOS_AXISYM_UPDATED_LAGRANGIAN_U_J_ELEMENT_H_INCLUDED
11 
12 // System includes
13 
14 // External includes
15 
16 // Project includes
18 
19 namespace Kratos
20 {
35 
36 
37 
38  class KRATOS_API(PFEM_SOLID_MECHANICS_APPLICATION) AxisymUpdatedLagrangianUJElement
40  {
41  public:
42 
48  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
53 
57 
60 
63 
65  AxisymUpdatedLagrangianUJElement(IndexType NewId, GeometryType::Pointer pGeometry);
66 
67  AxisymUpdatedLagrangianUJElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
68 
71 
72 
75 
79 
82 
83 
87 
98  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override ;
99 
107  Element::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override;
108 
109  //************* GETTING METHODS
110 
111 
112  //************* STARTING - ENDING METHODS
113 
118  void Initialize() override;
119 
120  //************* COMPUTING METHODS
121 
128  void CalculateMassMatrix(MatrixType& rMassMatrix,
129  ProcessInfo& rCurrentProcessInfo) override;
130 
131 
132  //************************************************************************************
133  //************************************************************************************
141  int Check(const ProcessInfo& rCurrentProcessInfo) override;
142 
146  virtual void CalculateKinematics(ElementDataType& rVariables,
147  const double& rPointNumber) override;
148 
149 
153  void CalculateDeformationGradient(const Matrix& rDN_DX,
154  Matrix& rF,
155  Matrix& rDeltaPosition,
156  double & rCurrentRadius,
157  double & rReferenceRadius);
158 
162  virtual void CalculateRadius(double & rCurrentRadius, double & rReferenceRadius, const Vector& rN);
163 
164 
168  void CalculateDeformationMatrix(Matrix& rB,
169  Matrix& rF,
170  Vector& rN,
171  double& rCurrentRadius);
172 
173 
177  void CalculateGreenLagrangeStrain(const Matrix& rF,
178  Vector& rStrainVector) override;
179 
183  void CalculateAlmansiStrain(const Matrix& rF,
184  Vector& rStrainVector) override;
185 
186 
187 
191 
202  protected:
208 
209 
213 
214 
218 
223  virtual void CalculateAndAddLHS(LocalSystemComponents& rLocalSystem,
224  ElementDataType& rVariables,
225  double& rIntegrationWeight) override;
226 
231  virtual void CalculateAndAddRHS(LocalSystemComponents& rLocalSystem,
232  ElementDataType& rVariables,
233  Vector& rVolumeForce,
234  double& rIntegrationWeight) override;
235 
239  virtual void InitializeElementData(ElementDataType & rVariables, const ProcessInfo& rCurrentProcessInfo) override;
240 
241 
242 
243 
247  virtual void CalculateAndAddKuum(MatrixType& rK,
248  ElementDataType & rVariables,
249  double& rIntegrationWeight
250  ) override;
251 
255  virtual void CalculateAndAddKuug(MatrixType& rK,
256  ElementDataType & rVariables,
257  double& rIntegrationWeight
258  ) override;
259 
263  virtual void CalculateAndAddKuJ (MatrixType& rK,
264  ElementDataType & rVariables,
265  double& rIntegrationWeight
266  ) override;
267 
271  virtual void CalculateAndAddKJu(MatrixType& rK,
272  ElementDataType & rVariables,
273  double& rIntegrationWeight
274  ) override;
275 
276 
277 
281  virtual void CalculateAndAddKJJ(MatrixType& rK,
282  ElementDataType & rVariables,
283  double& rIntegrationWeight
284  ) override;
285 
289  virtual void CalculateAndAddKJJStab(MatrixType& rK, ElementDataType & rVariables,
290  double& rIntegrationWeight
291  ) override;
292 
293 
297  virtual void CalculateAndAddJacobianForces(VectorType& rRightHandSideVector,
298  ElementDataType & rVariables,
299  double& rIntegrationWeight
300  ) override;
301 
305  virtual void CalculateAndAddStabilizedJacobian(VectorType& rRightHandSideVector,
306  ElementDataType & rVariables,
307  double& rIntegrationWeight
308  ) override;
309 
310 
311 
312 
316  virtual void GetHistoricalVariables( ElementDataType& rVariables,
317  const double& rPointNumber ) override;
318 
319 
320  /*
321  * Function to modify the deformation gradient to the constitutitve equation
322  */
323  virtual void ComputeConstitutiveVariables( ElementDataType& rVariables, Matrix& rFT, double& rdetFT) override;
324 
335 
336  private:
337 
343 
344 
348 
349 
353 
354 
359 
363  friend class Serializer;
364 
365  // A private default constructor necessary for serialization
366 
367  virtual void save(Serializer& rSerializer) const override;
368 
369  virtual void load(Serializer& rSerializer) override;
370 
371 
378 
379 
380 }; // Class AxisymUpdatedLagrangianUJElement
381 
382 
383 
384 } // namespace Kratos
385 #endif // KRATOS_UPDATED_LAGRANGIAN_U_J_ELEMENT_H_INCLUDED
386 
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Definition: axisym_updated_lagrangian_U_J_element.hpp:40
ConstitutiveLaw ConstitutiveLawType
Definition: axisym_updated_lagrangian_U_J_element.hpp:46
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(AxisymUpdatedLagrangianUJElement)
Counted pointer of LargeDisplacementUPElement.
ConstitutiveLawType::StressMeasure StressMeasureType
StressMeasure from constitutive laws.
Definition: axisym_updated_lagrangian_U_J_element.hpp:50
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: axisym_updated_lagrangian_U_J_element.hpp:52
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: axisym_updated_lagrangian_U_J_element.hpp:48
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
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
Updated Lagrangian Large Displacement Lagrangian U-wP Element for 3D and 2D geometries....
Definition: updated_lagrangian_U_J_element.hpp:44
double CalculateRadius(const PairedCondition *pCondition, const Vector &rNSlave)
Calculates the radius of axisymmetry.
Definition: mortar_explicit_contribution_utilities.cpp:675
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
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