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.
total_lagrangian.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: Riccardo Rossi
10 // Vicente Mataix Ferrandiz
11 //
12 
13 #pragma once
14 
15 // System includes
16 
17 // External include
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 {
34 
38 
42 
51 class KRATOS_API(STRUCTURAL_MECHANICS_APPLICATION) TotalLagrangian
52  : public BaseSolidElement
53 {
54 public:
57 
60 
62  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
63 
66 
69 
71  typedef std::size_t IndexType;
72 
74  typedef std::size_t SizeType;
75 
78 
82 
84  TotalLagrangian(IndexType NewId, GeometryType::Pointer pGeometry);
85  TotalLagrangian(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
86 
87  // Copy constructor
89  :BaseType(rOther)
90  {};
91 
93  ~TotalLagrangian() override;
94 
101 
109  Element::Pointer Create(
110  IndexType NewId,
111  GeometryType::Pointer pGeom,
112  PropertiesType::Pointer pProperties
113  ) const override;
114 
122  Element::Pointer Create(
123  IndexType NewId,
124  NodesArrayType const& ThisNodes,
125  PropertiesType::Pointer pProperties
126  ) const override;
127 
135  Element::Pointer Clone (
136  IndexType NewId,
137  NodesArrayType const& rThisNodes
138  ) const override;
139 
140  //std::string Info() const;
141 
142  void CalculateSensitivityMatrix(const Variable<array_1d<double, 3>>& rDesignVariable,
143  Matrix& rOutput,
144  const ProcessInfo& rCurrentProcessInfo) override;
145 
149 
156 
158  std::string Info() const override
159  {
160  std::stringstream buffer;
161  buffer << "Updated Lagrangian Solid Element #" << Id() << "\nConstitutive law: " << BaseType::mConstitutiveLawVector[0]->Info();
162  return buffer.str();
163  }
164 
166  void PrintInfo(std::ostream& rOStream) const override
167  {
168  rOStream << "Updated Lagrangian Solid Element #" << Id() << "\nConstitutive law: " << BaseType::mConstitutiveLawVector[0]->Info();
169  }
170 
172  void PrintData(std::ostream& rOStream) const override
173  {
174  pGetGeometry()->PrintData(rOStream);
175  }
176 
181 
182 protected:
188 
192 
194  {
195  }
196 
205  void CalculateAll(
206  MatrixType& rLeftHandSideMatrix,
207  VectorType& rRightHandSideVector,
208  const ProcessInfo& rCurrentProcessInfo,
209  const bool CalculateStiffnessMatrixFlag,
210  const bool CalculateResidualVectorFlag
211  ) override;
212 
219  void CalculateKinematicVariables(
220  KinematicVariables& rThisKinematicVariables,
221  const IndexType PointNumber,
222  const GeometryType::IntegrationMethod& rIntegrationMethod
223  ) override;
224 
228  std::size_t GetStrainSize() const;
229 
243 
244 private:
247 
251 
255 
259 
266  void CalculateB(Matrix& rB, Matrix const& rF, const Matrix& rDN_DX);
267 
268  void Calculate2DB(Matrix& rB, const Matrix& rF, const Matrix& rDN_DX);
269 
270  void Calculate3DB(Matrix& rB, const Matrix& rF, const Matrix& rDN_DX);
271 
272  void CalculateAxisymmetricB(Matrix& rB, const Matrix& rF, const Matrix& rDN_DX, const Vector& rN);
273 
274  void CalculateAxisymmetricF(Matrix const& rJ, Matrix const& rInvJ0, Vector const& rN, Matrix& rF);
275 
276  void CalculateStress(Vector& rStrain,
277  std::size_t IntegrationPoint,
278  Vector& rStress,
279  ProcessInfo const& rCurrentProcessInfo);
280 
281  void CalculateStress(Matrix const& rF,
282  std::size_t IntegrationPoint,
283  Vector& rStress,
284  ProcessInfo const& rCurrentProcessInfo);
285 
286  void CalculateStrain(Matrix const& rF,
287  std::size_t IntegrationPoint,
288  Vector& rStrain,
289  ProcessInfo const& rCurrentProcessInfo);
290 
291  void CalculateShapeSensitivity(ShapeParameter Deriv,
292  Matrix& rDN_DX0,
293  Matrix& rDN_DX0_Deriv,
294  Matrix& rF_Deriv,
295  double& rDetJ0_Deriv,
296  std::size_t IntegrationPointIndex);
297 
298  void CalculateBSensitivity(Matrix const& rDN_DX,
299  Matrix const& rF,
300  Matrix const& rDN_DX_Deriv,
301  Matrix const& rF_Deriv,
302  Matrix& rB_Deriv);
303 
304 
305  bool IsAxissymmetric() const;
306 
311 
315  friend class Serializer;
316 
317  // A private default constructor necessary for serialization
318 
319  void save(Serializer& rSerializer) const override;
320 
321  void load(Serializer& rSerializer) override;
322 
329  //TotalLagrangian& operator=(const TotalLagrangian& rOther);
331  //TotalLagrangian(const TotalLagrangian& rOther);
333 
334 }; // Class TotalLagrangian
335 
343 
344 } // 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 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
Short class definition.
Definition: integration_point.h:52
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
Total Lagrangian element for 2D and 3D geometries.
Definition: total_lagrangian.h:53
BaseSolidElement BaseType
The base element type.
Definition: total_lagrangian.h:68
std::size_t IndexType
The definition of the index type.
Definition: total_lagrangian.h:71
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(TotalLagrangian)
Counted pointer of TotalLagrangian.
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: total_lagrangian.h:172
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: total_lagrangian.h:62
std::size_t SizeType
The definition of the sizetype.
Definition: total_lagrangian.h:74
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: total_lagrangian.h:65
TotalLagrangian(TotalLagrangian const &rOther)
Definition: total_lagrangian.h:88
TotalLagrangian()
Definition: total_lagrangian.h:193
ConstitutiveLaw ConstitutiveLawType
Reference type definition for constitutive laws.
Definition: total_lagrangian.h:59
std::string Info() const override
Turn back information as a string.
Definition: total_lagrangian.h:158
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: total_lagrangian.h:166
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
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
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
def load(f)
Definition: ode_solve.py:307
Definition: geometrical_sensitivity_utility.h:33