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_segregated_V_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: April 2018 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_UPDATED_LAGRANGIAN_SEGREGATED_V_P_ELEMENT_H_INCLUDED)
11 #define KRATOS_UPDATED_LAGRANGIAN_SEGREGATED_V_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) UpdatedLagrangianSegregatedVPElement : public LargeDisplacementSegregatedVPElement
46 {
47 public:
48 
54  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
61 
63 
67 
70 
73 
75  UpdatedLagrangianSegregatedVPElement(IndexType NewId, GeometryType::Pointer pGeometry);
76 
77  UpdatedLagrangianSegregatedVPElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
78 
81 
82 
85 
89 
92 
96 
107  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
108 
116  Element::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override;
117 
118  //************* GETTING METHODS
119  //SET
123  void SetValuesOnIntegrationPoints(const Variable<double>& rVariable, const std::vector<double>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
124 
125  //GET:
129  void CalculateOnIntegrationPoints(const Variable<double>& rVariable, std::vector<double>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
130 
131  //************* STARTING - ENDING METHODS
136  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
137 
138  //************************************************************************************
139  //************************************************************************************
147  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
148 
152 
160  std::string Info() const override
161  {
162  std::stringstream buffer;
163  buffer << "Updated Lagrangian Segregated V-P Solid Element #" << Id();
164  return buffer.str();
165  }
166 
168  void PrintInfo(std::ostream& rOStream) const override
169  {
170  rOStream << "Updated Lagrangian Segregated V-P Solid Element #" << Id();
171  }
172 
174  void PrintData(std::ostream& rOStream) const override
175  {
176  GetGeometry().PrintData(rOStream);
177  }
182 
183 protected:
189 
193  std::vector< Matrix > mDeformationGradientF0;
194 
199 
206 
210  void FinalizeStepVariables(ElementDataType & rVariables,
211  const double& rPointNumber ) override;
212 
216  void InitializeElementData(ElementDataType & rVariables,
217  const ProcessInfo& rCurrentProcessInfo) override;
218 
222  void CalculateKinematics(ElementDataType& rVariables,
223  const double& rPointNumber) override;
224 
228  void CalculateKinetics(ElementDataType& rVariables,
229  const double& rPointNumber) override;
230 
234  void CalculateMaterialResponse(ElementDataType& rVariables,
236  const int & rPointNumber) override;
237 
241  void CalculateAndAddKpp(MatrixType& rLeftHandSideMatrix,
242  ElementDataType& rVariables) override;
243 
247  void CalculateAndAddPressureForces(VectorType& rRightHandSideVector,
248  ElementDataType & rVariables) override;
249 
253  void CalculateStabilizationTau(ElementDataType& rVariables);
254 
258  void GetFreeSurfaceFaces(std::vector<std::vector<SizeType> >& Faces);
259 
263  void GetFaceNormal(const std::vector<SizeType>& rFace, const ElementDataType & rVariables, Vector& rNormal);
264 
268  void GetFaceWeight(const std::vector<SizeType>& rFace, const ElementDataType & rVariables, double& rWeight, double& rNormalSize);
269 
273  void GetFaceNormal(const std::vector<SizeType>& rFace, Vector& rNormal);
274 
275 
279  void GetHistoricalVariables( ElementDataType& rVariables,
280  const double& rPointNumber ) override;
281 
285  double& CalculateVolumeChange(double& rVolumeChange, ElementDataType& rVariables) override;
286 
297 
298 private:
299 
318  friend class Serializer;
319 
320  // A private default constructor necessary for serialization
321 
322  void save(Serializer& rSerializer) const override;
323 
324  void load(Serializer& rSerializer) override;
325 
326 
333 
334 }; // Class UpdatedLagrangianSegregatedVPElement
335 
343 
344 } // namespace Kratos.
345 #endif // KRATOS_UPDATED_LAGRANGIAN_SEGREGATED_V_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 V Element for 3D and 2D geometries.
Definition: large_displacement_segregated_V_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
Updated Lagrangian Segregated Solid Element for 3D and 2D geometries.
Definition: updated_lagrangian_segregated_V_P_element.hpp:46
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: updated_lagrangian_segregated_V_P_element.hpp:54
GeometryData::SizeType SizeType
Type for size.
Definition: updated_lagrangian_segregated_V_P_element.hpp:60
std::vector< Matrix > mDeformationGradientF0
Definition: updated_lagrangian_segregated_V_P_element.hpp:193
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: updated_lagrangian_segregated_V_P_element.hpp:168
GlobalPointersVector< Element > ElementWeakPtrVectorType
Definition: updated_lagrangian_segregated_V_P_element.hpp:62
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: updated_lagrangian_segregated_V_P_element.hpp:58
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(UpdatedLagrangianSegregatedVPElement)
Counted pointer of UpdatedLagrangianSegregatedVPElement.
ConstitutiveLaw ConstitutiveLawType
Definition: updated_lagrangian_segregated_V_P_element.hpp:52
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: updated_lagrangian_segregated_V_P_element.hpp:174
Vector mDeterminantF0
Definition: updated_lagrangian_segregated_V_P_element.hpp:198
ConstitutiveLawType::StressMeasure StressMeasureType
StressMeasure from constitutive laws.
Definition: updated_lagrangian_segregated_V_P_element.hpp:56
std::string Info() const override
Turn back information as a string.
Definition: updated_lagrangian_segregated_V_P_element.hpp:160
void CalculateStabilizationTau(double &rTau, double &rElementLength, const array_1d< double, TDim > &rVelocity, const Matrix &rContravariantMetricTensor, const double Reaction, const double EffectiveKinematicViscosity, const double Alpha, const double Gamma, const double DeltaTime, const double DynamicTau)
Definition: convection_diffusion_reaction_stabilization_utilities.h:52
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: constitutive_law.h:189
Definition: solid_element.hpp:83