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.
linear_solid_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_LINEAR_SOLID_ELEMENT_H_INCLUDED)
11 #define KRATOS_LINEAR_SOLID_ELEMENT_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
18 #include "includes/checks.h"
19 #include "includes/element.h"
20 #include "custom_utilities/solid_mechanics_math_utilities.hpp"
21 
22 namespace Kratos
23 {
38 
40 
41 class KRATOS_API(SOLID_MECHANICS_APPLICATION) LinearSolidElement
42  : public Element
43 {
44 public:
45 
51  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
58 
62 
63 
64 public:
65 
66 
69 
71  LinearSolidElement(IndexType NewId, GeometryType::Pointer pGeometry);
72 
73  LinearSolidElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
74 
77 
79  ~LinearSolidElement() override;
80 
84 
87 
91 
102  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
103 
111  Element::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override;
112 
113  //************* GETTING METHODS
114 
119  IntegrationMethod GetIntegrationMethod() const override;
120 
124  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
125 
129  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
130 
134  void GetValuesVector(Vector& rValues, int Step = 0) const override;
135 
139  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
140 
144  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
145 
146 
147 
148  //on integration points:
149  // see element.h
150 
154  void CalculateOnIntegrationPoints(const Variable<Matrix>& rVariable, std::vector<Matrix>& rValues, const ProcessInfo& rCurrentProcessInfo) override;
155 
156 
157  //************* STARTING - ENDING METHODS
158 
163  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
164 
168  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
169 
173  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
174 
178  void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
179 
183  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
184 
185 
186  //************* COMPUTING METHODS
187 
197  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix,
198  VectorType& rRightHandSideVector,
199  const ProcessInfo& rCurrentProcessInfo) override;
200 
201 
208  void CalculateMassMatrix(MatrixType& rMassMatrix,
209  const ProcessInfo& rCurrentProcessInfo) override;
210 
217  void CalculateDampingMatrix(MatrixType& rDampingMatrix,
218  const ProcessInfo& rCurrentProcessInfo) override;
219 
220 
230  void AddExplicitContribution(const VectorType& rRHSVector,
231  const Variable<VectorType>& rRHSVariable,
232  const Variable<array_1d<double,3> >& rDestinationVariable,
233  const ProcessInfo& rCurrentProcessInfo) override;
234 
235 
236 
237  //************************************************************************************
238  //************************************************************************************
246  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
247 
248 
262 
263 protected:
269 
274 
278  std::vector<ConstitutiveLaw::Pointer> mConstitutiveLawVector;
279 
280 
285 
286  //constexpr const std::size_t& Dimension() const {return GetGeometry().WorkingSpaceDimension();}
287 
291 
292 
296  void ClearNodalForces ();
297 
298 
302  Matrix& CalculateTotalDeltaPosition(Matrix & rDeltaPosition);
303 
304 
308  void CalculateInfinitesimalStrain(Vector& rStrainVector,
309  const Matrix& rDN_DX);
310 
314  Vector& CalculateVolumeForce(Vector& rVolumeForce, const Vector& rN);
315 
316 
327 
328 private:
329 
348  friend class Serializer;
349 
350  // A private default constructor necessary for serialization
351 
352  void save(Serializer& rSerializer) const override;
353 
354  void load(Serializer& rSerializer) override;
355 
356 
363 
364 }; // Class LinearSolidElement
365 
373 
374 } // namespace Kratos.
375 #endif // KRATOS_LINEAR_SOLID_ELEMENT_H_INCLUDED defined
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
Definition: constitutive_law.h:47
StressMeasure
Definition: constitutive_law.h:69
Base class for all Elements.
Definition: element.h:60
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
std::size_t IndexType
Definition: flags.h:74
IntegrationMethod
Definition: geometry_data.h:76
std::size_t SizeType
Definition: geometry_data.h:173
Linear Solid Element for 3D and 2D geometries. (Template for learning element design)
Definition: linear_solid_element.hpp:43
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: linear_solid_element.hpp:55
ConstitutiveLaw ConstitutiveLawType
Definition: linear_solid_element.hpp:49
IntegrationMethod mThisIntegrationMethod
Definition: linear_solid_element.hpp:273
std::vector< ConstitutiveLaw::Pointer > mConstitutiveLawVector
Definition: linear_solid_element.hpp:278
GeometryData::SizeType SizeType
Type for size.
Definition: linear_solid_element.hpp:57
ConstitutiveLawType::StressMeasure StressMeasureType
StressMeasure from constitutive laws.
Definition: linear_solid_element.hpp:53
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(LinearSolidElement)
Counted pointer of LinearSolidElement.
LinearSolidElement()
Definition: linear_solid_element.hpp:284
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: linear_solid_element.hpp:51
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
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
void InitializeSolutionStep(ConstructionUtility &rThisUtil, std::string ThermalSubModelPartName, std::string MechanicalSubModelPartName, std::string HeatFluxSubModelPartName, std::string HydraulicPressureSubModelPartName, bool thermal_conditions, bool mechanical_conditions, int phase)
Definition: add_custom_utilities_to_python.cpp:45
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
def load(f)
Definition: ode_solve.py:307