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.
large_displacement_beam_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: November 2015 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_LARGE_DISPLACEMENT_BEAM_ELEMENT_H_INCLUDED)
11 #define KRATOS_LARGE_DISPLACEMENT_BEAM_ELEMENT_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
19 
20 namespace Kratos
21 {
36 
38 
46 class KRATOS_API(SOLID_MECHANICS_APPLICATION) LargeDisplacementBeamElement
47  :public BeamElement
48 {
49 public:
50 
56  typedef ConstitutiveLawType::Pointer ConstitutiveLawPointerType;
69 
72 
76 
78  LargeDisplacementBeamElement(IndexType NewId, GeometryType::Pointer pGeometry);
79 
80  LargeDisplacementBeamElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
81 
84 
86  ~LargeDisplacementBeamElement() override;
87 
88 
92 
93 
101  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
102 
103 
108  void Initialize(const ProcessInfo& rCurrentProcessInfo) override;
109 
113  void InitializeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
114 
118  void InitializeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
119 
123  void FinalizeNonLinearIteration(const ProcessInfo& rCurrentProcessInfo) override;
124 
125 
129  void FinalizeSolutionStep(const ProcessInfo& rCurrentProcessInfo) override;
130 
131 
132  //************* COMPUTING METHODS
133 
134 
141  void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) override;
142 
143 
144  //on integration points:
148  void CalculateOnIntegrationPoints(const Variable<double>& rVariable,
149  std::vector<double>& rOutput,
150  const ProcessInfo& rCurrentProcessInfo) override;
151 
152 
157  std::vector< array_1d<double, 3 > >& Output,
158  const ProcessInfo& rCurrentProcessInfo) override;
159 
160 
161  //************************************************************************************
162  //************************************************************************************
170  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
171 
181 
183  std::string Info() const override
184  {
185  std::stringstream buffer;
186  buffer << "Large Displacement Beam Element #" << Id();
187  return buffer.str();
188  }
189 
191  void PrintInfo(std::ostream& rOStream) const override
192  {
193  rOStream << "Large Displacement Beam Element #" << Id();
194  }
195 
197  void PrintData(std::ostream& rOStream) const override
198  {
199  GetGeometry().PrintData(rOStream);
200  }
201 
206 
207 protected:
208 
214 
219 
224 
229 
233  double mInvJ0;
234 
238  std::vector<Vector> mCurrentCurvatureVectors;
239 
243  std::vector<Vector> mPreviousCurvatureVectors;
244 
248  std::vector<QuaternionType> mFrameQuaternionsReduced;
249 
253  std::vector<QuaternionType> mFrameQuaternionsFull;
254 
255 
260 
264 
269  IntegrationMethod GetReducedIntegrationMethod() const;
270 
275  IntegrationMethod GetFullIntegrationMethod() const;
276 
280  void CalculateDynamicSystem( LocalSystemComponents& rLocalSystem,
281  const ProcessInfo& rCurrentProcessInfo ) override;
282 
286  virtual void MapToSpatialFrame(const ElementDataType& rVariables, Matrix& rVariable);
287 
288 
292  virtual void CalculateCurrentCurvature(ElementDataType& rVariables, const Variable<array_1d<double, 3 > >& rVariable);
293 
294 
298  void CalculateKinematics(ElementDataType& rVariables,
299  const unsigned int& rPointNumber) override;
300 
304  virtual void CalculateFrameMapping(ElementDataType& rVariables,
305  const unsigned int& rPointNumber);
306 
307 
311  virtual void UpdateStrainVariables(ElementDataType& rVariables,
312  const unsigned int& rPointNumber);
313 
314 
318  void CalculateConstitutiveMatrix(ElementDataType& rVariables) override;
319 
320 
324  void CalculateStressResultants(ElementDataType& rVariables, const unsigned int& rPointNumber) override;
325 
330  void CalculateAndAddKuum(MatrixType& rLeftHandSideMatrix,
331  ElementDataType& rVariables,
332  double& rIntegrationWeight) override;
333 
334 
335 
339  void CalculateAndAddKuug(MatrixType& rLeftHandSideMatrix,
340  ElementDataType& rVariables,
341  double& rIntegrationWeight) override;
342 
343  virtual void CalculateAndAddKuug2(MatrixType& rLeftHandSideMatrix,
344  ElementDataType& rVariables,
345  double& rIntegrationWeight);
346 
350  virtual void CalculateAndAddKuuf(MatrixType& rLeftHandSideMatrix,
351  ElementDataType& rVariables,
352  double& rIntegrationWeight);
353 
354 
355 
359  virtual void CalculateAndAddFollowerForces(VectorType& rRightHandSideVector,
360  ElementDataType& rVariables,
361  double& rIntegrationWeight);
362 
363 
367  void CalculateAndAddExternalForces(VectorType& rRightHandSideVector,
368  ElementDataType& rVariables,
369  Vector& rVolumeForce,
370  double& rIntegrationWeight) override;
371 
372 
373 
377  void CalculateAndAddInertiaLHS(MatrixType& rLeftHandSideMatrix,
378  ElementDataType& rVariables,
379  const ProcessInfo& rCurrentProcessInfo,
380  double& rIntegrationWeight) override;
381 
385  void CalculateAndAddInertiaRHS(VectorType& rRightHandSideVector,
386  ElementDataType& rVariables,
387  const ProcessInfo& rCurrentProcessInfo,
388  double& rIntegrationWeight) override;
389 
393  void CalculateAndAddInternalForces(VectorType& rRightHandSideVector,
394  ElementDataType & rVariables,
395  double& rIntegrationWeight) override;
396 
397 
401  void CalculateOperator(MatrixType& rOperator,
402  Vector& rN,
403  const int& rNode);
404 
408  void CalculateDiscreteOperator(MatrixType& rDiscreteOperator,
409  ElementDataType& rVariables,
410  const int& rNode);
411 
415  virtual void CalculateDifferentialOperator(MatrixType& rDifferentialOperator,
416  ElementDataType& rVariables,
417  const int& rNode,
418  double alpha = 0);
419 
420 
424  void CalculateBmatrix(MatrixType& rBmatrix,
425  ElementDataType& rVariables);
426 
430  double& CalculateTotalMass( SectionProperties& Section, double& rTotalMass );
431 
435  void CalculateInertiaDyadic(SectionProperties& rSection, Matrix& rInertiaDyadic);
436 
440  virtual void CalculateRotationLinearPartTensor(Vector& rRotationVector, Matrix& rRotationTensor);
441 
442 
446  void GetCurrentNodalMovements(Vector& rValues, const int& rNode, unsigned int PointNumber = 0);
447 
451  void GetCurrentNodalVelocities(Vector& rValues, const int& rNode, unsigned int PointNumber = 0);
452 
456  void GetKineticMatrix(Matrix& rKineticMatrix, const double& rMass, const Matrix& rInertia);
457 
458 
462  virtual void CalculateExternalForcesEnergy(double& rEnergy, ElementDataType& rVariables, Vector& rVolumeForce, double& rIntegrationWeight);
463 
467  virtual void CalculateInternalForcesEnergy(double& rEnergy, ElementDataType& rVariables, double& rIntegrationWeight);
468 
469 
473  virtual void CalculateStrainEnergy(double& rEnergy, ElementDataType& rVariables, const ProcessInfo& rCurrentProcessInfo, double& rIntegrationWeight);
474 
478  virtual void CalculateKineticEnergy(double& rEnergy, ElementDataType& rVariables, const ProcessInfo& rCurrentProcessInfo, double& rIntegrationWeight);
479 
483  virtual void CalculateMomentumRelations(array_1d<double,3>& LinearMomentum, array_1d<double,3>& AngularMomentum, ElementDataType& rVariables, const ProcessInfo& rCurrentProcessInfo, double& rIntegrationWeight);
484 
495 
496 private:
497 
516  friend class Serializer;
517 
518 
519  // A private default constructor necessary for serialization
520 
521 
522  void save(Serializer& rSerializer) const override;
523 
524  void load(Serializer& rSerializer) override;
525 
532 
533 
534 }; // Class LargeDisplacementBeamElement
535 
536 } // namespace Kratos.
537 #endif // KRATOS_LARGE_DISPLACEMENT_BEAM_ELEMENT_H_INCLUDED defined
Beam Element for 2D and 3D space dimensions.
Definition: beam_element.hpp:44
Definition: beam_math_utilities.hpp:31
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
Beam Element for 3D space dimension Simo Displacement-Rotation Geometrically Exact Rod element.
Definition: large_displacement_beam_element.hpp:48
IntegrationMethod mReducedIntegrationMethod
Definition: large_displacement_beam_element.hpp:218
ConstitutiveLawType::StressMeasure StressMeasureType
StressMeasure from constitutive laws.
Definition: large_displacement_beam_element.hpp:58
std::vector< Vector > mCurrentCurvatureVectors
Definition: large_displacement_beam_element.hpp:238
GeometryData::IntegrationMethod IntegrationMethod
Type definition for integration methods.
Definition: large_displacement_beam_element.hpp:60
LargeDisplacementBeamElement()
Definition: large_displacement_beam_element.hpp:259
std::vector< QuaternionType > mFrameQuaternionsReduced
Definition: large_displacement_beam_element.hpp:248
double mInvJ0
Definition: large_displacement_beam_element.hpp:233
ConstitutiveLaw ConstitutiveLawType
Definition: large_displacement_beam_element.hpp:54
IntegrationMethod mFullIntegrationMethod
Definition: large_displacement_beam_element.hpp:223
std::string Info() const override
Turn back information as a string.
Definition: large_displacement_beam_element.hpp:183
Quaternion< double > QuaternionType
Type definition for quaternion.
Definition: large_displacement_beam_element.hpp:64
int mIterationCounter
Definition: large_displacement_beam_element.hpp:228
BeamMathUtils< double > BeamMathUtilsType
Type definition for beam utilities.
Definition: large_displacement_beam_element.hpp:62
std::vector< Vector > mPreviousCurvatureVectors
Definition: large_displacement_beam_element.hpp:243
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: large_displacement_beam_element.hpp:197
BeamElement::ElementDataType ElementDataType
Type for element variables.
Definition: large_displacement_beam_element.hpp:68
GeometryData::SizeType SizeType
Type for size.
Definition: large_displacement_beam_element.hpp:66
std::vector< QuaternionType > mFrameQuaternionsFull
Definition: large_displacement_beam_element.hpp:253
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(LargeDisplacementBeamElement)
Counted pointer of LargeDisplacementBeamElement.
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: large_displacement_beam_element.hpp:191
ConstitutiveLawType::Pointer ConstitutiveLawPointerType
Pointer type for constitutive laws.
Definition: large_displacement_beam_element.hpp:56
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
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
double CalculateKineticEnergy(Element &rElement)
Definition: mpm_energy_calculation_utility.cpp:60
double CalculateStrainEnergy(Element &rElement)
Definition: mpm_energy_calculation_utility.cpp:89
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
alpha
Definition: generate_convection_diffusion_explicit_element.py:113
def load(f)
Definition: ode_solve.py:307
Definition: beam_element.hpp:123