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_Pressure_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_UPDATED_LAGRANGIAN_UP_SECOND_ELEMENT_H_INCLUDED )
11 #define KRATOS_UPDATED_LAGRANGIAN_UP_SECOND_ELEMENT_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
19 
20 namespace Kratos
21 {
36 
38 
39  // THIS IS THE SECOND VERSION OF THE UP-ELEMENT FOR ANY CONSTITUTIVE EQUATION.
40  // In this version, the pressure is the pressure and it is not some "2D" invented pressure,..
41  // the constitutive equation should not be the UP version
42 
43 
44  class KRATOS_API(PFEM_SOLID_MECHANICS_APPLICATION) UpdatedLagrangianUPressureElement
46  {
47 
48  protected:
49  typedef struct
50  {
51  unsigned int voigtsize;
52 
55 
57 
59 
61 
62  public:
63 
69  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
78 
82 
85 
88 
90  UpdatedLagrangianUPressureElement(IndexType NewId, GeometryType::Pointer pGeometry);
91 
92  UpdatedLagrangianUPressureElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
93 
96 
97 
100 
104 
107 
108 
112 
123  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
124 
132  Element::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override;
133 
134 
135 
136 
137  //************************************************************************************
138  //************************************************************************************
146  int Check(const ProcessInfo& rCurrentProcessInfo) override;
147 
148  /*
149  * Get on rVariable a Matrix Value from the Element Constitutive Law
150  */
151  void GetValueOnIntegrationPoints(const Variable<Matrix>& rVariable, std::vector<Matrix>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
152 
153 
154  // Get Values defined in order to avoid a clang warning (?)
155  void GetValueOnIntegrationPoints(const Variable<double>& rVariable, std::vector<double>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
156 
160  virtual void CalculateOnIntegrationPoints(const Variable<Matrix >& rVariable, std::vector< Matrix >& rOutput, const ProcessInfo& rCurrentProcessInfo) override;
161 
162 
166 
170 
174 
178 
180  protected:
183 
187 
189 
191 
195 
199 
204  virtual void CalculateAndAddLHS(LocalSystemComponents& rLocalSystem,
205  ElementDataType& rVariables,
206  double& rIntegrationWeight) override;
207 
212  virtual void CalculateAndAddRHS(LocalSystemComponents& rLocalSystem,
213  ElementDataType& rVariables,
214  Vector& rVolumeForce,
215  double& rIntegrationWeight) override;
216 
220  virtual void InitializeElementData(ElementDataType & rVariables,
221  const ProcessInfo& rCurrentProcessInfo) override;
222 
223 
224 
228  virtual void CalculateAndAddKuumElemUP(MatrixType& rK,
229  ElementDataType & rVariables,
230  ThisElementData & rElementVariables,
231  double& rIntegrationWeight
232  );
233 
237  virtual void CalculateAndAddKuugElemUP(MatrixType& rK,
238  ElementDataType & rVariables,
239  ThisElementData & rElementVariables,
240  double& rIntegrationWeight
241  );
242 
246  virtual void CalculateAndAddKupElemUP (MatrixType& rK,
247  ElementDataType & rVariables,
248  ThisElementData & rElementVariables,
249  double& rIntegrationWeight
250  );
251 
255  virtual void CalculateAndAddKpuElemUP(MatrixType& rK,
256  ElementDataType & rVariables,
257  ThisElementData & rElementVariables,
258  double& rIntegrationWeight
259  );
260 
261 
265  virtual void CalculateAndAddKppElemUP(MatrixType& rK,
266  ElementDataType & rVariables,
267  ThisElementData & rElementVariables,
268  double& rIntegrationWeight
269  );
270 
271 
275  virtual void CalculateAndAddKppStabElemUP(MatrixType& rK,
276  ElementDataType & rVariables,
277  ThisElementData & rElementVariables,
278  double& rIntegrationWeight
279  );
280 
281  /*
282  * Calculation of the Internal Forces due to sigma. Fi = B * sigma
283  */
284  void CalculateAndAddInternalForcesElemUP(VectorType& rRightHandSideVector,
285  ElementDataType & rVariables,
286  ThisElementData & rElementVariables,
287  double& rIntegrationWeight
288  );
289 
293  virtual void CalculateAndAddPressureForcesElemUP(VectorType& rRightHandSideVector,
294  ElementDataType & rVariables,
295  ThisElementData & rElementVariables,
296  double& rIntegrationWeight
297  );
298 
299 
303  virtual void CalculateAndAddStabilizedPressureElemUP(VectorType& rRightHandSideVector,
304  ElementDataType & rVariables,
305  ThisElementData & rElementVariables,
306  double& rIntegrationWeight
307  );
308 
312  virtual void CalculateThisElementData( ThisElementData& rElementVariables, const ElementDataType & rVariables);
313 
314 
325 
326  private:
327 
333 
334 
338 
339 
343 
344 
349 
353  friend class Serializer;
354 
355  // A private default constructor necessary for serialization
356 
357  virtual void save(Serializer& rSerializer) const override;
358 
359  virtual void load(Serializer& rSerializer) override;
360 
361 
368 
369 
370 }; // Class UpdatedLagrangianUPressureElement
371 
372 
373 
374 } // namespace Kratos
375 #endif // KRATOS_____
376 
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
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
Large Displacement Lagrangian U-P Element for 3D and 2D geometries. Linear Triangles and Tetrahedra (...
Definition: updated_lagrangian_U_Pressure_element.hpp:46
ConstitutiveLawType::StressMeasure StressMeasureType
StressMeasure from constitutive laws.
Definition: updated_lagrangian_U_Pressure_element.hpp:71
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: updated_lagrangian_U_Pressure_element.hpp:69
double mElementStabilizationNumber
Definition: updated_lagrangian_U_Pressure_element.hpp:188
GeometryData::SizeType SizeType
Type for size.
Definition: updated_lagrangian_U_Pressure_element.hpp:75
ConstitutiveLaw ConstitutiveLawType
Definition: updated_lagrangian_U_Pressure_element.hpp:67
double mElementScalingNumber
Definition: updated_lagrangian_U_Pressure_element.hpp:190
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: updated_lagrangian_U_Pressure_element.hpp:73
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(UpdatedLagrangianUPressureElement)
Counted pointer of LargeDisplacementUPElement.
UpdatedLagrangianUPElement::ElementDataType ElementDataType
Type for element variables.
Definition: updated_lagrangian_U_Pressure_element.hpp:77
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
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
Definition: updated_lagrangian_U_Pressure_element.hpp:50
unsigned int voigtsize
Definition: updated_lagrangian_U_Pressure_element.hpp:51
Matrix DeviatoricTensor
Definition: updated_lagrangian_U_Pressure_element.hpp:58
double ElementalMeanStress
Definition: updated_lagrangian_U_Pressure_element.hpp:54
double NodalMeanStress
Definition: updated_lagrangian_U_Pressure_element.hpp:53
Vector StressVector
Definition: updated_lagrangian_U_Pressure_element.hpp:56