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.
small_displacement_element.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosDamApplication $
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_SMALL_DISPLACEMENT_ELEMENT_H_INCLUDED)
11 #define KRATOS_SMALL_DISPLACEMENT_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(DAM_APPLICATION) SmallDisplacementElement
46  : public SolidElement
47 {
48 public:
49 
52 
56  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
65 
68 
70 
71 
74 
77 
79  SmallDisplacementElement(IndexType NewId, GeometryType::Pointer pGeometry);
80 
81  SmallDisplacementElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
82 
85 
87  ~SmallDisplacementElement() override;
88 
92 
95 
99 
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 
119  //on integration points:
120 
124  void CalculateOnIntegrationPoints(const Variable<double>& rVariable,
125  std::vector<double>& rOutput,
126  const ProcessInfo& rCurrentProcessInfo) override;
127 
131  void CalculateOnIntegrationPoints(const Variable<Vector>& rVariable,
132  std::vector<Vector>& rOutput,
133  const ProcessInfo& rCurrentProcessInfo) override;
134 
138  void CalculateOnIntegrationPoints(const Variable<Matrix >& rVariable,
139  std::vector< Matrix >& rOutput,
140  const ProcessInfo& rCurrentProcessInfo) override;
141 
142 
143  //************************************************************************************
144  //************************************************************************************
152  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
153 
157 
168 
169 protected:
172 
176 
180 
184 
189  void CalculateAndAddKuug(MatrixType& rLeftHandSideMatrix,
190  ElementDataType& rVariables,
191  double& rIntegrationWeight) override;
192 
193 
197  void SetElementData(ElementDataType& rVariables,
199  const int & rPointNumber) override;
200 
204  void CalculateKinematics(ElementDataType& rVariables,
205  const double& rPointNumber) override;
206 
210  void InitializeElementData(ElementDataType & rVariables,
211  const ProcessInfo& rCurrentProcessInfo) override;
212 
213 
217  virtual void CalculateInfinitesimalStrain(const Matrix& rH,
218  Vector& rStrainVector);
219 
223  void CalculateDisplacementGradient(Matrix& rH,
224  const Matrix& rDN_DX);
225 
226 
237 
238 private:
239 
245 
246 
250 
251 
255 
256 
261 
265  friend class Serializer;
266 
267  // A private default constructor necessary for serialization
268 
269  void save(Serializer& rSerializer) const override;
270 
271  void load(Serializer& rSerializer) override;
272 
273 
280 
281 }; // Class SmallDisplacementElement
282 
290 
291 } // namespace Kratos.
292 #endif // KRATOS_SMALL_DISPLACEMENT_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
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
Small Displacement Element for 3D and 2D geometries.
Definition: small_displacement_element.hpp:47
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: small_displacement_element.hpp:60
GeometryData::SizeType SizeType
Type for size.
Definition: small_displacement_element.hpp:62
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: small_displacement_element.hpp:56
ConstitutiveLaw ConstitutiveLawType
Reference type definition for constitutive laws.
Definition: small_displacement_element.hpp:54
ConstitutiveLawType::StressMeasure StressMeasureType
StressMeasure from constitutive laws.
Definition: small_displacement_element.hpp:58
SolidElement::ElementDataType ElementDataType
Type for element variables.
Definition: small_displacement_element.hpp:64
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(SmallDisplacementElement)
Counted pointer of SmallDisplacementElement.
Large Displacement Lagrangian Element for 3D and 2D geometries. (base class)
Definition: solid_element.hpp:49
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: constitutive_law.h:189
Definition: solid_element.hpp:83