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.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: August 2017 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_SMALL_DISPLACEMENT_BEAM_ELEMENT_H_INCLUDED)
11 #define KRATOS_SMALL_DISPLACEMENT_BEAM_ELEMENT_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
19 
20 namespace Kratos
21 {
36 
38 
43 class KRATOS_API(SOLID_MECHANICS_APPLICATION) SmallDisplacementBeamElement
44  :public BeamElement
45 {
46 public:
47 
50 
53 
57 
59  SmallDisplacementBeamElement(IndexType NewId, GeometryType::Pointer pGeometry);
60 
61  SmallDisplacementBeamElement(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties);
62 
65 
67  ~SmallDisplacementBeamElement() override;
68 
69 
73 
81  Element::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override;
82 
83 
84  //************************************************************************************
85  //************************************************************************************
93  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
97 
105  std::string Info() const override
106  {
107  std::stringstream buffer;
108  buffer << "Small Displacement Beam Element #" << Id();
109  return buffer.str();
110  }
111 
113  void PrintInfo(std::ostream& rOStream) const override
114  {
115  rOStream << "Small Displacement Beam Element #" << Id();
116  }
117 
119  void PrintData(std::ostream& rOStream) const override
120  {
121  GetGeometry().PrintData(rOStream);
122  }
127 
128 protected:
129 
139 
140  //constexpr const std::size_t& Dimension() const {return GetGeometry().WorkingSpaceDimension();}
141 
145 
149  void InitializeElementData(ElementDataType & rVariables, const ProcessInfo& rCurrentProcessInfo) override;
150 
154  void CalculateKinematics(ElementDataType& rVariables,
155  const unsigned int& rPointNumber) override;
156 
160  virtual void CalculateDeformationMatrix(Matrix& rB, const Vector& rN, const Matrix& rDN_DX);
161 
162 
166  void CalculateConstitutiveMatrix(ElementDataType& rVariables) override;
167 
168 
172  void CalculateStressResultants(ElementDataType& rVariables, const unsigned int& rPointNumber) override;
173 
177  void CalculateTransformationMatrix(Matrix& RotationMatrix);
178 
179 
183  void MapLocalToGlobal(ElementDataType& rVariables, Matrix& rVariable) override;
184 
188  void MapLocalToGlobal(ElementDataType& rVariables, VectorType& rVector) override;
189 
190 
191 
196  void CalculateAndAddKuum(MatrixType& rLeftHandSideMatrix,
197  ElementDataType& rVariables,
198  double& rIntegrationWeight) override;
199 
200 
204  void CalculateLocalStiffnessMatrix(Matrix& LocalMatrix,ElementDataType& rVariables);
205 
206 
210  void CalculateAndAddInternalForces(VectorType& rRightHandSideVector,
211  ElementDataType & rVariables,
212  double& rIntegrationWeight) override;
213 
214 
225 
226 private:
227 
246  friend class Serializer;
247 
248 
249  // A private default constructor necessary for serialization
250 
251 
252  void save(Serializer& rSerializer) const override
253  {
255  }
256 
257  void load(Serializer& rSerializer) override
258  {
260  }
261 
268 
269 
270 }; // Class SmallDisplacementBeamElement
271 
272 } // namespace Kratos.
273 #endif // KRATOS_SMALL_DISPLACEMENT_BEAM_ELEMENT_H_INCLUDED defined
Beam Element for 2D and 3D space dimensions.
Definition: beam_element.hpp:44
Base class for all Elements.
Definition: element.h:60
std::size_t IndexType
Definition: flags.h:74
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 displacements beam element for 2D and 3D space.
Definition: small_displacement_beam_element.hpp:45
SmallDisplacementBeamElement()
Definition: small_displacement_beam_element.hpp:138
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: small_displacement_beam_element.hpp:119
std::string Info() const override
Turn back information as a string.
Definition: small_displacement_beam_element.hpp:105
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(SmallDisplacementBeamElement)
Counted pointer of SmallDisplacementBeamElement.
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: small_displacement_beam_element.hpp:113
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
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
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
def load(f)
Definition: ode_solve.py:307