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.
swimming_particle.h
Go to the documentation of this file.
1 //
2 // Project Name: Kratos
3 // Last Modified by: $Author: G.Casas, gcasas@cimne.upc.edu $
4 // Date: $Date: 2019-1-02 16:07:33 $
5 // Revision: $Revision: 1.1.1.1 $
6 //
7 //
8 #if !defined(KRATOS_SWIMMING_PARTICLE_H_INCLUDED)
9 #define KRATOS_SWIMMING_PARTICLE_H_INCLUDED
10 
11 // System includes
12 #include <string>
13 #include <iostream>
14 
15 // Project includes
16 #include "includes/define.h"
17 #include "includes/model_part.h"
22 
23 namespace Kratos
24 {
25 
26  template< class TBaseElement >
27  class KRATOS_API(SWIMMING_DEM_APPLICATION) SwimmingParticle : public TBaseElement
28  {
29  public:
30 
31  typedef std::size_t IndexType;
32  typedef Node NodeType;
36 
37  using TBaseElement::GetGeometry;
39  using TBaseElement::mRealMass;
40  using TBaseElement::mRadius;
41  using TBaseElement::CalculateVolume;
42  using TBaseElement::GetMass;
43  using TBaseElement::GetForce;
44 
45 
48 
51 
55 
57  SwimmingParticle():TBaseElement(){}
58  SwimmingParticle(IndexType NewId, GeometryType::Pointer pGeometry):TBaseElement(NewId, pGeometry){}
59  SwimmingParticle(IndexType NewId, NodesArrayType const& ThisNodes):TBaseElement(NewId, ThisNodes){}
60  SwimmingParticle(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
61  :TBaseElement(NewId, pGeometry, pProperties){}
62 
63  Element::Pointer Create(IndexType NewId,
64  NodesArrayType const& ThisNodes,
65  PropertiesType::Pointer pProperties) const override
66  {
67  return Element::Pointer(new SwimmingParticle(NewId, GetGeometry().Create(ThisNodes), pProperties));
68  };
69 
71  virtual ~SwimmingParticle(){}
72 
74  virtual void CreateHydrodynamicInteractionLaws(const ProcessInfo& r_process_info);
75 
76  void ComputeAdditionalForces(array_1d<double, 3>& additionally_applied_force,
77  array_1d<double, 3>& additionally_applied_moment,
78  const ProcessInfo& rCurrentProcessInfo,
79  const array_1d<double,3>& gravity) override;
80 
81  std::vector<Node::Pointer> mNeighbourNodes;
82  std::vector<double> mNeighbourNodesDistances;
83 
85  virtual std::string Info() const override
86  {
87  std::stringstream buffer;
88  buffer << "Swimming version of " << TBaseElement::Info();
89  return buffer.str();
90  }
91 
93  virtual void PrintInfo(std::ostream& rOStream) const override
94  {
95  rOStream << "Swimming version of " << TBaseElement::Info();
96  }
97 
99  virtual void PrintData(std::ostream& rOStream) const override{}
100 
101  void Calculate(const Variable<array_1d<double, 3 > >& rVariable,
102  array_1d<double, 3 > & Output,
103  const ProcessInfo& r_current_process_info) override;
104 
105  protected:
106 
107  void MemberDeclarationFirstStep(const ProcessInfo& r_current_process_info) override;
108  void AdditionalCalculate(const Variable<double>& rVariable, double& Output, const ProcessInfo& r_current_process_info) override;
109  array_1d<double,3> ComputeWeight(const array_1d<double,3>& gravity, const ProcessInfo& r_process_info) override;
110  void AddCentrifugalForces(array_1d<double,3>& weight, const ProcessInfo& r_process_info);
111  void AddCoriolisForces(array_1d<double,3>& weight, const ProcessInfo& r_process_info);
112  void AddRelativeAccelerationForces(array_1d<double,3>& weight, const ProcessInfo& r_process_info);
113  void AddEulerForces(array_1d<double,3>& weight, const ProcessInfo& r_process_info);
114  virtual double GetFluidMass();
115 
118 
119 
123 
124 
128 
129 
133 
134 
138 
139 
143 
144 
148 
149 
151 
152  private:
153  double CalculateDragCoeffFromSphericity(const double reynolds, double sphericity, int drag_modifier_type);
154  void UpdateNodalValues(NodeType& node,
155  const array_1d<double, 3>& hydrodynamic_force,
156  const array_1d<double, 3>& hydrodynamic_moment,
157  const array_1d<double, 3>& weight,
158  const array_1d<double, 3>& buoyancy,
159  const array_1d<double, 3>& drag_force,
160  const array_1d<double, 3>& inviscid_force,
161  const array_1d<double, 3>& history_force,
162  const array_1d<double, 3>& vorticity_induced_lift,
163  const array_1d<double, 3>& rotation_induced_lift,
164  const double &force_reduction_coeff,
165  const ProcessInfo& r_current_process_info);
166  void ApplyNumericalAveragingWithOldForces(NodeType& node,
167  array_1d<double, 3>& additionally_applied_force,
168  const ProcessInfo& r_current_process_info);
169  void ApplyDragPorosityModification(double& drag_coeff);
170  void Initialize(const ProcessInfo& r_process_info) override;
173 
174 
178 
179  bool mFirstStep;
180  int mPorosityCorrectionType;
181  double mFluidDensity;
182  double mKinematicViscosity;
183  double mSphericity;
184  double mNormOfSlipVel;
185  array_1d<double, 3> mSlipVel;
186  HydrodynamicInteractionLaw::Pointer mHydrodynamicInteractionLaw;
187 
191 
192 
196 
197 
201 
202 
206 
207 
211 
212 
216 
217  friend class Serializer;
218 
219  virtual void save(Serializer& rSerializer) const override
220  {
222  }
223 
224  virtual void load(Serializer& rSerializer) override
225  {
227  }
228 
229 
230 
232 
234 
235 
237 
238  }; // Class SwimmingParticle
239 
241 
244 
245 
249 
250 } // namespace Kratos.
251 
252 #endif // KRATOS_SWIMMING_PARTICLE_H_INCLUDED defined
253 
254 
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
std::string Info() const override
Turn back information as a string.
Definition: periodic_interface_process.hpp:93
Definition: discrete_element.h:38
Geometry base class.
Definition: geometry.h:71
This class defines the node.
Definition: node.h:65
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
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Definition: swimming_particle.h:28
virtual void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: swimming_particle.h:93
Element::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Definition: swimming_particle.h:63
SwimmingParticle()
Default constructor.
Definition: swimming_particle.h:57
virtual void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: swimming_particle.h:99
SwimmingParticle(IndexType NewId, NodesArrayType const &ThisNodes)
Definition: swimming_particle.h:59
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(SwimmingParticle)
Pointer definition of SwimmingParticle.
std::vector< Node::Pointer > mNeighbourNodes
Definition: swimming_particle.h:81
Node NodeType
Definition: swimming_particle.h:32
virtual ~SwimmingParticle()
Destructor.
Definition: swimming_particle.h:71
std::vector< double > mNeighbourNodesDistances
Definition: swimming_particle.h:82
std::size_t IndexType
Definition: swimming_particle.h:31
Geometry< NodeType > GeometryType
Definition: swimming_particle.h:33
Geometry< NodeType >::PointsArrayType NodesArrayType
Definition: swimming_particle.h:34
Properties PropertiesType
Definition: swimming_particle.h:35
SwimmingParticle(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Definition: swimming_particle.h:60
SwimmingParticle(IndexType NewId, GeometryType::Pointer pGeometry)
Definition: swimming_particle.h:58
virtual std::string Info() const override
Turn back information as a string.
Definition: swimming_particle.h:85
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
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
TDataType Calculate(GeometryType &dummy, const Variable< TDataType > &rVariable)
Definition: add_geometries_to_python.cpp:103
double GetDensity(const Properties &rProps, const IndexType Index)
Definition: shell_utilities.cpp:223
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
def load(f)
Definition: ode_solve.py:307
Definition: mesh_converter.cpp:38