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_beam_element_3D2N.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosSolidMechanicsApplication $
3 // Created by: $Author: Nelson Maireni $
4 // Last modified by: $Co-Author: JMCarbonell $
5 // Date: $Date: November 2015 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_SMALL_DISPLACEMENT_BEAM_ELEMENT_3D2N_H_INCLUDED)
11 #define KRATOS_SMALL_DISPLACEMENT_BEAM_ELEMENT_3D2N_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 
21 namespace Kratos
22 {
37 
39 
45 class KRATOS_API(SOLID_MECHANICS_APPLICATION) SmallDisplacementBeamElement3D2N
46  :public Element
47 {
48 public:
49 
55 
58 
60 
61 protected:
62 
66  KRATOS_DEFINE_LOCAL_FLAG( COMPUTE_RHS_VECTOR );
67  KRATOS_DEFINE_LOCAL_FLAG( COMPUTE_LHS_MATRIX );
68 
70  {
71 
72  double Area; // Area or the beam section
73  double Inertia_z; // Moment of Inertia about the local z axis, Iz local
74  double Inertia_y; // Moment of Inertia about the local y axis, Iy local
75  double Polar_Inertia; // Polar Moment of Inertia, measures an object's ability to resist twisting, when acted upon by differences of torque along its length.
76  };
77 
86  {
87  private:
88 
89  //for calculation local system with compacted LHS and RHS
90  MatrixType *mpLeftHandSideMatrix;
91  VectorType *mpRightHandSideVector;
92 
93  public:
94 
95  //calculation flags
98 
99 
103  void SetLeftHandSideMatrix( MatrixType& rLeftHandSideMatrix ) { mpLeftHandSideMatrix = &rLeftHandSideMatrix; };
104 
105  void SetRightHandSideVector( VectorType& rRightHandSideVector ) { mpRightHandSideVector = &rRightHandSideVector; };
106 
107 
111  MatrixType& GetLeftHandSideMatrix() { return *mpLeftHandSideMatrix; };
112 
113  VectorType& GetRightHandSideVector() { return *mpRightHandSideVector; };
114 
115  };
116 
117 
118 public:
119 
122 
124  SmallDisplacementBeamElement3D2N(IndexType NewId, GeometryType::Pointer pGeometry);
125 
126  SmallDisplacementBeamElement3D2N(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
127 
130 
133 
134 
138 
139 
147  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
148 
149 
150 
151 
152  //************* GETTING METHODS
153 
159  IntegrationMethod GetIntegrationMethod() const override;
160 
164  void GetDofList(DofsVectorType& rElementalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
165 
169  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
170 
174  void GetValuesVector(Vector& rValues, int Step = 0) const override;
175 
179  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
180 
184  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
185 
186  //on integration points:
198  //************* STARTING - ENDING METHODS
199 
204  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
205 
209  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
210 
211 
215  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
216 
217  //************* COMPUTING METHODS
218 
219 
228  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
229 
236  void CalculateRightHandSide(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
237 
238 
245  void CalculateLeftHandSide(MatrixType& rLeftHandSideMatrix, const ProcessInfo& rCurrentProcessInfo) override;
246 
247 
254  void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) override;
255 
256 
257  //on integration points:
262  std::vector< array_1d<double, 3 > >& Output,
263  const ProcessInfo& rCurrentProcessInfo) override;
264 
265 
266  //************************************************************************************
267  //************************************************************************************
275  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
276 
280 
288  std::string Info() const override
289  {
290  std::stringstream buffer;
291  buffer << "Beam Element #" << Id();
292  return buffer.str();
293  }
294 
296  void PrintInfo(std::ostream& rOStream) const override
297  {
298  rOStream << "Beam Element #" << Id();
299  }
300 
302  void PrintData(std::ostream& rOStream) const override
303  {
304  GetGeometry().PrintData(rOStream);
305  }
310 
311 protected:
312 
318 
323 
324 
329 
330 
334  double mLength; // Length of the beam element
335 
340 
341  //constexpr const std::size_t& Dimension() const {return GetGeometry().WorkingSpaceDimension();}
342 
346 
347 
354  void CalculateElementalSystem( LocalSystemComponents& rLocalSystem,
355  const ProcessInfo& rCurrentProcessInfo );
356 
357 
361  void CalculateSectionProperties();
362 
366  virtual void InitializeSystemMatrices(MatrixType& rLeftHandSideMatrix,
367  VectorType& rRightHandSideVector,
368  Flags& rCalculationFlags);
369 
373  void CalculateAndAddLHS(LocalSystemComponents& rLocalSystem);
374 
378  void CalculateAndAddRHS(LocalSystemComponents& rLocalSystem, Vector& VolumeForce);
379 
383  void CalculateLocalStiffnessMatrix(Matrix& LocalMatrix);
384 
385 
389  void CalculateTransformationMatrix(Matrix& RotationMatrix);
390 
391 
395  virtual Vector& CalculateVolumeForce(Vector& rVolumeForce, const Vector& rN);
396 
401  void CalculateGlobalBodyForce(Vector& rGlobalForceVector, Vector& rVolumeForce);
402 
403 
404  void CalculateLocalBodyForce(Vector& rLocalBody, Vector& rVolumeForce);
405 
406 
407  void CalculateDistributedBodyForce(const int Direction, Vector& Load, Vector& rVolumeForce);
408 
409 
413  double CalculateInternalAxil ( const double& N, const double& AxialLoad, const double& X );
414 
415  double CalculateInternalShear ( const double& Q, const double& ShearLoad, const double& X );
416 
417  double CalculateInternalMoment( const double& M1, const double& M2, const double& ShearLoad, const double& X );
418 
419 
423  void CalculateLocalNodalStress(Vector& Stress, Vector& rVolumeForce );
424 
435 
436 private:
437 
456  friend class Serializer;
457 
458 
459  // A private default constructor necessary for serialization
460 
461 
462  void save(Serializer& rSerializer) const override
463  {
465  }
466 
467  void load(Serializer& rSerializer) override
468  {
470  }
471 
478 
479 
480 }; // Class SmallDisplacementBeamElement3D2N
481 
482 } // namespace Kratos.
483 #endif // KRATOS_SMALL_DISPLACEMENT_BEAM_ELEMENT_3D2N_H_INCLUDED defined
Base class for all Elements.
Definition: element.h:60
Definition: flags.h:58
IntegrationMethod
Definition: geometry_data.h:76
std::size_t SizeType
Definition: geometry_data.h:173
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
Beam Element for 3D space dimension.
Definition: small_displacement_beam_element_3D2N.hpp:47
std::string Info() const override
Turn back information as a string.
Definition: small_displacement_beam_element_3D2N.hpp:288
GeometryData::IntegrationMethod IntegrationMethod
Definition: small_displacement_beam_element_3D2N.hpp:52
GeometryData::SizeType SizeType
Type for size.
Definition: small_displacement_beam_element_3D2N.hpp:54
SmallDisplacementBeamElement3D2N()
Definition: small_displacement_beam_element_3D2N.hpp:339
KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_LHS_MATRIX)
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: small_displacement_beam_element_3D2N.hpp:302
SectionProperties mSection
Definition: small_displacement_beam_element_3D2N.hpp:328
IntegrationMethod mThisIntegrationMethod
Definition: small_displacement_beam_element_3D2N.hpp:322
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: small_displacement_beam_element_3D2N.hpp:296
KRATOS_DEFINE_LOCAL_FLAG(COMPUTE_RHS_VECTOR)
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(SmallDisplacementBeamElement3D2N)
Counted pointer of SmallDisplacementBeamElement3D2N.
double mLength
Definition: small_displacement_beam_element_3D2N.hpp:334
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
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 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
ModelPart::NodesContainerType NodesArrayType
Definition: gid_gauss_point_container.h:42
Stress
Definition: isotropic_damage_automatic_differentiation.py:135
tuple Q
Definition: isotropic_damage_automatic_differentiation.py:235
def load(f)
Definition: ode_solve.py:307
N
Definition: sensitivityMatrix.py:29
Definition: small_displacement_beam_element_3D2N.hpp:86
VectorType & GetRightHandSideVector()
Definition: small_displacement_beam_element_3D2N.hpp:113
void SetLeftHandSideMatrix(MatrixType &rLeftHandSideMatrix)
Definition: small_displacement_beam_element_3D2N.hpp:103
void SetRightHandSideVector(VectorType &rRightHandSideVector)
Definition: small_displacement_beam_element_3D2N.hpp:105
MatrixType RotationMatrix
Definition: small_displacement_beam_element_3D2N.hpp:97
Flags CalculationFlags
Definition: small_displacement_beam_element_3D2N.hpp:96
MatrixType & GetLeftHandSideMatrix()
Definition: small_displacement_beam_element_3D2N.hpp:111
Definition: small_displacement_beam_element_3D2N.hpp:70
double Inertia_y
Definition: small_displacement_beam_element_3D2N.hpp:74
double Area
Definition: small_displacement_beam_element_3D2N.hpp:72
double Polar_Inertia
Definition: small_displacement_beam_element_3D2N.hpp:75
double Inertia_z
Definition: small_displacement_beam_element_3D2N.hpp:73