3 #if !defined KRATOS_RIGID_BODY_ELEMENT_H_INCLUDED
4 #define KRATOS_RIGID_BODY_ELEMENT_H_INCLUDED
48 void Initialize(
const ProcessInfo& r_process_info)
override;
49 virtual void SetIntegrationScheme(DEMIntegrationScheme::Pointer& translational_integration_scheme, DEMIntegrationScheme::Pointer& rotational_integration_scheme);
50 virtual void CustomInitialize(
ModelPart& rigid_body_element_sub_model_part);
52 virtual void UpdateLinearDisplacementAndVelocityOfNodes();
53 virtual void UpdateAngularDisplacementAndVelocityOfNodes();
55 virtual void CollectForcesAndTorquesFromTheNodes();
59 virtual void Move(
const double delta_t,
const bool rotation_option,
const double force_reduction_factor,
const int StepFlag);
63 virtual double GetMass();
67 virtual std::string
Info()
const override
69 std::stringstream buffer;
70 buffer <<
"Discrete Element #" << Id();
75 virtual void PrintInfo(std::ostream& rOStream)
const override
77 rOStream <<
"Discrete Element #" << Id();
81 virtual void PrintData(std::ostream& rOStream)
const override
100 virtual void save(
Serializer& rSerializer)
const override;
113 rOStream << std::endl;
Definition: dem_integration_scheme.h:24
Base class for all Elements.
Definition: element.h:60
virtual void Initialize(const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:341
std::size_t IndexType
Definition: flags.h:74
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
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
Definition: rigid_body_element.h:31
virtual std::string Info() const override
Turn back information as a string.
Definition: rigid_body_element.h:67
DEMIntegrationScheme * mpTranslationalIntegrationScheme
Definition: rigid_body_element.h:88
virtual DEMIntegrationScheme & GetRotationalIntegrationScheme()
Definition: rigid_body_element.h:61
virtual DEMIntegrationScheme & GetTranslationalIntegrationScheme()
Definition: rigid_body_element.h:60
virtual void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: rigid_body_element.h:75
std::vector< array_1d< double, 3 > > mListOfCoordinates
Definition: rigid_body_element.h:86
std::vector< Node::Pointer > mListOfNodes
Definition: rigid_body_element.h:87
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(RigidBodyElement3D)
Pointer definition of RigidBodyElement3D.
std::vector< RigidFace3D * > mListOfRigidFaces
Definition: rigid_body_element.h:93
array_1d< double, 3 > mInertias
Definition: rigid_body_element.h:94
double GetSqrtOfRealMass()
virtual void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: rigid_body_element.h:81
DEMIntegrationScheme * mpRotationalIntegrationScheme
Definition: rigid_body_element.h:89
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
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
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
float velocity
Definition: PecletTest.py:54
def load(f)
Definition: ode_solve.py:307
float delta_t
Definition: rotatingcone_PureConvectionBenchmarking.py:129