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.
cr_beam_element_linear_3D2N.hpp
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: Klaus B. Sautter
10 //
11 
12 #pragma once
13 
14 // System includes
15 
16 // External includes
17 
18 // Project includes
20 #include "includes/define.h"
21 #include "includes/variables.h"
22 
23 namespace Kratos
24 {
33 class KRATOS_API(STRUCTURAL_MECHANICS_APPLICATION) CrBeamElementLinear3D2N : public CrBeamElement3D2N
34 {
35 
36 public:
38 
40  CrBeamElementLinear3D2N(IndexType NewId, GeometryType::Pointer pGeometry);
41  CrBeamElementLinear3D2N(IndexType NewId, GeometryType::Pointer pGeometry,
42  PropertiesType::Pointer pProperties);
43 
44 
45  ~CrBeamElementLinear3D2N() override;
46 
54  Element::Pointer Create(
55  IndexType NewId,
56  GeometryType::Pointer pGeom,
57  PropertiesType::Pointer pProperties
58  ) const override;
59 
67  Element::Pointer Create(
68  IndexType NewId,
69  NodesArrayType const& ThisNodes,
70  PropertiesType::Pointer pProperties
71  ) const override;
72 
73  void CalculateLocalSystem(
74  MatrixType& rLeftHandSideMatrix,
75  VectorType& rRightHandSideVector,
76  const ProcessInfo& rCurrentProcessInfo) override;
77 
78  void CalculateRightHandSide(
79  VectorType& rRightHandSideVector,
80  const ProcessInfo& rCurrentProcessInfo) override;
81 
82  void CalculateLeftHandSide(
83  MatrixType& rLeftHandSideMatrix,
84  const ProcessInfo& rCurrentProcessInfo) override;
85 
86  void CalculateMassMatrix(
87  MatrixType& rMassMatrix,
88  const ProcessInfo& rCurrentProcessInfo) override;
89 
93  BoundedMatrix<double,msLocalSize,msLocalSize> CalculateDeformationStiffness() const override;
94 
95  void Calculate(const Variable<Matrix>& rVariable, Matrix& rOutput,
96  const ProcessInfo& rCurrentProcessInfo) override;
97 
99  const Variable<array_1d<double, 3 > >& rVariable,
100  std::vector< array_1d<double, 3 > >& rOutput,
101  const ProcessInfo& rCurrentProcessInfo) override;
102 
104  const Variable<Vector >& rVariable,
105  std::vector< Vector >& rOutput,
106  const ProcessInfo& rCurrentProcessInfo) override;
107 
108 private:
109 
110  friend class Serializer;
111  void save(Serializer& rSerializer) const override;
112  void load(Serializer& rSerializer) override;
113 };
114 
115 }
This is a 3D-2node beam element with 3 translational dofs and 3 rotational dof per node.
Definition: cr_beam_element_3D2N.hpp:35
This is a linear 3D-2node beam element with 3 translational dofs and 3 rotational dof per node inheri...
Definition: cr_beam_element_linear_3D2N.hpp:34
CrBeamElementLinear3D2N()
Definition: cr_beam_element_linear_3D2N.hpp:39
void CalculateOnIntegrationPoints(const Variable< array_1d< double, 3 > > &rVariable, std::vector< array_1d< double, 3 > > &rOutput, const ProcessInfo &rCurrentProcessInfo) override
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(CrBeamElementLinear3D2N)
Definition: amatrix_interface.h:41
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
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
TDataType Calculate(GeometryType &dummy, const Variable< TDataType > &rVariable)
Definition: add_geometries_to_python.cpp:103
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