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_bbar.h
Go to the documentation of this file.
1 // KRATOS ___| | | |
2 // \___ \ __| __| | | __| __| | | __| _` | |
3 // | | | | | ( | | | | ( | |
4 // _____/ \__|_| \__,_|\___|\__|\__,_|_| \__,_|_| MECHANICS
5 //
6 // License: BSD License
7 // license: StructuralMechanicsApplication/license.txt
8 //
9 // Main authors: Marcelo Raschi
10 // Manuel Caicedo
11 // Javier Mroginski
12 
13 #pragma once
14 
15 // System includes
16 
17 // External includes
18 
19 // Project includes
20 #include "includes/define.h"
21 #include "custom_elements/base_solid_element.h"
22 #include "includes/variables.h"
23 
24 namespace Kratos
25 {
40 
55 class KRATOS_API(STRUCTURAL_MECHANICS_APPLICATION) SmallDisplacementBbar
56  : public BaseSolidElement
57 {
58 protected:
63  : public KinematicVariables
64 {
66 
74  const SizeType StrainSize,
75  const SizeType Dimension,
76  const SizeType NumberOfNodes
77  ) : KinematicVariables(StrainSize, Dimension, NumberOfNodes)
78  {
79  Bh = ZeroVector(Dimension * NumberOfNodes);
80  }
81 };
82 
83 public:
89  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
92 
95 
97  typedef std::size_t IndexType;
98 
100  typedef std::size_t SizeType;
101 
104 
108 
110  SmallDisplacementBbar(IndexType NewId, GeometryType::Pointer pGeometry);
111  SmallDisplacementBbar(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
112 
113  // Copy constructor
115  :BaseType(rOther)
116  {};
117 
118 
120  ~SmallDisplacementBbar() override;
121 
128 
136  Element::Pointer Create(
137  IndexType NewId,
138  GeometryType::Pointer pGeom,
139  PropertiesType::Pointer pProperties
140  ) const override;
141 
149  Element::Pointer Create(
150  IndexType NewId,
151  NodesArrayType const& ThisNodes,
152  PropertiesType::Pointer pProperties
153  ) const override;
154 
162  const Variable<Matrix>& rVariable,
163  std::vector<Matrix>& rOutput,
164  const ProcessInfo& rCurrentProcessInfo
165  ) override;
166 
174  const Variable<Vector>& rVariable,
175  std::vector<Vector>& rOutput,
176  const ProcessInfo& rCurrentProcessInfo
177  ) override;
178 
186  const Variable<double>& rVariable,
187  std::vector<double>& rOutput,
188  const ProcessInfo& rCurrentProcessInfo
189  ) override;
190 
195  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
196 
197 protected:
198 
204 
208 
210  {
211  }
212 
216  bool UseElementProvidedStrain() const override;
217 
223  void CalculateKinematicVariables(
224  KinematicVariables& rThisKinematicVariables,
225  const IndexType PointNumber,
226  const GeometryType::IntegrationMethod& rIntegrationMethod
227  ) override;
228 
232  void CalculateAndAddResidualVector(
233  VectorType& rRightHandSideVector,
234  const KinematicVariables& rThisKinematicVariables,
235  const ProcessInfo& rCurrentProcessInfo,
236  const array_1d<double, 3>& rBodyForce,
237  const Vector& rStressVector,
238  const double IntegrationWeight
239  ) const override;
240 
250  void SetConstitutiveVariables(
251  KinematicVariables& rThisKinematicVariables,
252  ConstitutiveVariables& rThisConstitutiveVariables,
254  const IndexType PointNumber,
255  const GeometryType::IntegrationPointsArrayType& IntegrationPoints
256  ) override;
257 
258 
267  void CalculateAll(
268  MatrixType& rLeftHandSideMatrix,
269  VectorType& rRightHandSideVector,
270  const ProcessInfo& rCurrentProcessInfo,
271  const bool CalculateStiffnessMatrixFlag,
272  const bool CalculateResidualVectorFlag
273  ) override;
274 
280  void CalculateB(
281  Matrix& rB,
282  const Matrix& DN_DX
283  );
284 
285  Matrix ComputeEquivalentF(const Vector& rStrainTensor);
286 
292  void CalculateKinematicVariablesBbar(
293  KinematicVariablesBbar& rThisKinematicVariables,
294  const IndexType PointNumber,
295  const GeometryType::IntegrationPointsArrayType& IntegrationPoints
296  );
297 
303  void CalculateBbar(
304  Matrix &rB,
305  Vector &rBh,
306  const Matrix &DN_DX,
307  const GeometryType::IntegrationPointsArrayType &IntegrationPoints,
308  const IndexType PointNumber
309  );
310 
311  // Compute Bbar components
316  void CalculateHydrostaticDeformationMatrix(KinematicVariablesBbar& rThisKinematicVariables);
317 
318 private:
319 
320  friend class Serializer;
321 
322  // A private default constructor necessary for serialization
323 
324  void save(Serializer& rSerializer) const override;
325 
326  void load(Serializer& rSerializer) override;
327 
334  //SmallDisplacementStrElement& operator=(const SmallDisplacementStrElement& rOther);
336  //SmallDisplacementStrElement(const SmallDisplacementStrElement& rOther);
338 
339 }; // Class SmallDisplacementBbar
340 
348 
349 } // namespace Kratos.
This is base class used to define the solid elements.
Definition: base_solid_element.h:67
Definition: constitutive_law.h:47
std::size_t SizeType
Definition: element.h:94
std::size_t IndexType
Definition: flags.h:74
This defines the geometrical object, base definition of the element and condition entities.
Definition: geometrical_object.h:58
IntegrationMethod
Definition: geometry_data.h:76
std::vector< IntegrationPointType > IntegrationPointsArrayType
Definition: geometry.h:161
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
Infinitesimal strain definition with mixed B-bar formulation.
Definition: small_displacement_bbar.h:57
SmallDisplacementBbar()
Definition: small_displacement_bbar.h:209
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: small_displacement_bbar.h:89
std::size_t SizeType
The definition of the sizetype.
Definition: small_displacement_bbar.h:100
ConstitutiveLaw ConstitutiveLawType
Definition: small_displacement_bbar.h:87
BaseSolidElement BaseType
The base element type.
Definition: small_displacement_bbar.h:94
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: small_displacement_bbar.h:91
SmallDisplacementBbar(SmallDisplacementBbar const &rOther)
Definition: small_displacement_bbar.h:114
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(SmallDisplacementBbar)
Counted pointer of SmallDisplacementStrElement.
std::size_t IndexType
The definition of the index type.
Definition: small_displacement_bbar.h:97
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
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 ComputeEquivalentF(const Element &rElement, const TVectorType &rStrainTensor, TMatrixType &rF)
This method computes the deformation gradient F (for small deformation solid elements)
Definition: structural_mechanics_element_utilities.h:69
void CalculateB(const GeometricalObject &rElement, const TMatrixType1 &rDN_DX, TMatrixType2 &rB)
This method computes the deformation tensor B (for small deformation solid elements)
Definition: structural_mechanics_element_utilities.h:105
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
def load(f)
Definition: ode_solve.py:307
Definition: base_solid_element.h:73
Definition: constitutive_law.h:189
Definition: small_displacement_bbar.h:64
KinematicVariablesBbar(const SizeType StrainSize, const SizeType Dimension, const SizeType NumberOfNodes)
Definition: small_displacement_bbar.h:73
Vector Bh
Definition: small_displacement_bbar.h:65