SPHinXsys  alpha version
base_particle_dynamics.h
Go to the documentation of this file.
1 /* -------------------------------------------------------------------------*
2 * SPHinXsys *
3 * --------------------------------------------------------------------------*
4 * SPHinXsys (pronunciation: s'finksis) is an acronym from Smoothed Particle *
5 * Hydrodynamics for industrial compleX systems. It provides C++ APIs for *
6 * physical accurate simulation and aims to model coupled industrial dynamic *
7 * systems including fluid, solid, multi-body dynamics and beyond with SPH *
8 * (smoothed particle hydrodynamics), a meshless computational method using *
9 * particle discretization. *
10 * *
11 * SPHinXsys is partially funded by German Research Foundation *
12 * (Deutsche Forschungsgemeinschaft) DFG HU1527/6-1, HU1527/10-1 *
13 * and HU1527/12-1. *
14 * *
15 * Portions copyright (c) 2017-2020 Technical University of Munich and *
16 * the authors' affiliations. *
17 * *
18 * Licensed under the Apache License, Version 2.0 (the "License"); you may *
19 * not use this file except in compliance with the License. You may obtain a *
20 * copy of the License at http://www.apache.org/licenses/LICENSE-2.0. *
21 * *
22 * --------------------------------------------------------------------------*/
31 #ifndef BASE_PARTICLE_DYNAMICS_H
32 #define BASE_PARTICLE_DYNAMICS_H
33 
34 #include "base_data_package.h"
35 #include "sph_data_containers.h"
36 #include "neighbor_relation.h"
37 #include "all_body_relations.h"
38 #include "base_body.h"
39 
40 #include <functional>
41 
42 using namespace std::placeholders;
43 
44 namespace SPH
45 {
46 
47 
49  typedef std::function<void(size_t, Real)> ParticleFunctor;
51  template <class ReturnType>
52  using ReduceFunctor = std::function<ReturnType(size_t, Real)>;
53 
55  void ParticleIterator(size_t total_real_particles, const ParticleFunctor &particle_functor, Real dt = 0.0);
57  void ParticleIterator_parallel(size_t total_real_particles, const ParticleFunctor &particle_functor, Real dt = 0.0);
58 
60  template <class ReturnType, typename ReduceOperation>
61  ReturnType ReduceIterator(size_t total_real_particles, ReturnType temp,
62  ReduceFunctor<ReturnType> &reduce_functor, ReduceOperation &reduce_operation, Real dt = 0.0);
64  template <class ReturnType, typename ReduceOperation>
65  ReturnType ReduceIterator_parallel(size_t total_real_particles, ReturnType temp,
66  ReduceFunctor<ReturnType> &reduce_functor, ReduceOperation &reduce_operation, Real dt = 0.0);
67 
69  void ParticleIteratorSplittingSweep(SplitCellLists &split_cell_lists,
70  const ParticleFunctor &particle_functor, Real dt = 0.0);
73  const ParticleFunctor &particle_functor, Real dt = 0.0);
74 
76  template <class ReturnType>
77  struct ReduceSum
78  {
79  ReturnType operator()(const ReturnType &x, const ReturnType &y) const { return x + y; };
80  };
82  struct ReduceMax
83  {
84  Real operator()(Real x, Real y) const { return SMAX(x, y); };
85  };
87  struct ReduceMin
88  {
89  Real operator()(Real x, Real y) const { return SMIN(x, y); };
90  };
92  struct ReduceOR
93  {
94  bool operator()(bool x, bool y) const { return x || y; };
95  };
97  struct ReduceAND
98  {
99  bool operator()(bool x, bool y) const { return x && y; };
100  };
103  {
104  Vecd operator()(const Vecd &x, const Vecd &y) const
105  {
106  Vecd lower_bound;
107  for (int i = 0; i < lower_bound.size(); ++i)
108  lower_bound[i] = SMIN(x[i], y[i]);
109  return lower_bound;
110  };
111  };
114  {
115  Vecd operator()(const Vecd &x, const Vecd &y) const
116  {
117  Vecd upper_bound;
118  for (int i = 0; i < upper_bound.size(); ++i)
119  upper_bound[i] = SMAX(x[i], y[i]);
120  return upper_bound;
121  };
122  };
123 
129  {
130  public:
131  explicit GlobalStaticVariables(){};
132  virtual ~GlobalStaticVariables(){};
133 
135  static Real physical_time_;
136  };
137 
144  template <class ReturnType = void>
146  {
147  public:
148  explicit ParticleDynamics(SPHBody &sph_body);
149  virtual ~ParticleDynamics(){};
150 
151  SPHBody *getSPHBody() { return sph_body_; };
154  virtual ReturnType exec(Real dt = 0.0) = 0;
155  virtual ReturnType parallel_exec(Real dt = 0.0) = 0;
156 
157  protected:
158  SPHBody *sph_body_;
159  SPHAdaptation *sph_adaptation_;
160  BaseParticles *base_particles_;
161 
162  void setBodyUpdated() { sph_body_->setNewlyUpdated(); };
164  virtual void setupDynamics(Real dt = 0.0){};
165  };
166 
172  {
173  public:
174  explicit DataDelegateEmptyBase(SPHBody &sph_body){};
175  virtual ~DataDelegateEmptyBase(){};
176  };
177 
182  template <class BodyType = SPHBody,
183  class ParticlesType = BaseParticles,
184  class MaterialType = BaseMaterial>
186  {
187  public:
188  explicit DataDelegateSimple(SPHBody &sph_body)
189  : body_(DynamicCast<BodyType>(this, &sph_body)),
190  particles_(DynamicCast<ParticlesType>(this, sph_body.base_particles_)),
191  material_(DynamicCast<MaterialType>(this, sph_body.base_material_)),
192  sorted_id_(sph_body.base_particles_->sorted_id_),
193  unsorted_id_(sph_body.base_particles_->unsorted_id_){};
194  virtual ~DataDelegateSimple(){};
195 
196  BodyType *getBody() { return body_; };
197  ParticlesType *getParticles() { return particles_; };
198  MaterialType *getMaterial() { return material_; };
199 
200  protected:
201  BodyType *body_;
202  ParticlesType *particles_;
203  MaterialType *material_;
204  StdLargeVec<size_t> &sorted_id_;
205  StdLargeVec<size_t> &unsorted_id_;
206  };
207 
212  template <class BodyType = SPHBody,
213  class ParticlesType = BaseParticles,
214  class MaterialType = BaseMaterial,
217  {
218  public:
219  explicit DataDelegateInner(BaseBodyRelationInner &body_inner_relation)
220  : BaseDataDelegateType(*body_inner_relation.sph_body_),
221  inner_configuration_(body_inner_relation.inner_configuration_){};
222  virtual ~DataDelegateInner(){};
223 
224  protected:
226  ParticleConfiguration &inner_configuration_;
227  };
228 
233  template <class BodyType = SPHBody,
234  class ParticlesType = BaseParticles,
235  class MaterialType = BaseMaterial,
236  class ContactBodyType = SPHBody,
237  class ContactParticlesType = BaseParticles,
238  class ContactMaterialType = BaseMaterial,
239  class BaseDataDelegateType = DataDelegateSimple<BodyType, ParticlesType, MaterialType>>
241  {
242  public:
243  explicit DataDelegateContact(BaseBodyRelationContact &body_contact_relation);
244  virtual ~DataDelegateContact(){};
245 
246  protected:
247  StdVec<ContactBodyType *> contact_bodies_;
248  StdVec<ContactParticlesType *> contact_particles_;
249  StdVec<ContactMaterialType *> contact_material_;
251  StdVec<ParticleConfiguration *> contact_configuration_;
252  };
253 
258  template <class BodyType = SPHBody,
259  class ParticlesType = BaseParticles,
260  class MaterialType = BaseMaterial,
261  class ContactBodyType = SPHBody,
262  class ContactParticlesType = BaseParticles,
263  class ContactMaterialType = BaseMaterial>
264  class DataDelegateComplex : public DataDelegateInner<BodyType, ParticlesType, MaterialType>,
265  public DataDelegateContact<BodyType, ParticlesType, MaterialType,
266  ContactBodyType, ContactParticlesType, ContactMaterialType, DataDelegateEmptyBase>
267  {
268  public:
269  explicit DataDelegateComplex(ComplexBodyRelation &body_complex_relation)
270  : DataDelegateInner<BodyType, ParticlesType, MaterialType>(body_complex_relation.inner_relation_),
271  DataDelegateContact<BodyType, ParticlesType, MaterialType, ContactBodyType, ContactParticlesType,
272  ContactMaterialType, DataDelegateEmptyBase>(body_complex_relation.contact_relation_){};
273  virtual ~DataDelegateComplex(){};
274  };
275 
280  template <class ParticleDynamicsInnerType, class ContactDataType>
282  {
283  public:
285  BaseBodyRelationContact &contact_relation)
286  : ParticleDynamicsInnerType(inner_relation), ContactDataType(contact_relation){};
287 
289  BaseBodyRelationContact &extra_contact_relation);
290 
291  virtual ~ParticleDynamicsComplex(){};
292 
293  protected:
294  virtual void prepareContactData() = 0;
295  };
296 }
297 #endif //BASE_PARTICLE_DYNAMICS_H
StdVec< ConcurrentCellLists > SplitCellLists
Definition: sph_data_containers.h:54
Base class for all adaptations. The base class defines essential global parameters. It is also used for single-resolution method. In the constructor parameter, system_refinement_ratio defines the relation between present resolution to the system reference resolution. The derived classes are defined for more complex adaptations.
Definition: adaptation.h:55
Particles with essential (geometric and kinematic) data. There are three types of particles, all par...
Definition: base_particles.h:81
Definition: base_particle_dynamics.h:87
Definition: base_particle_dynamics.h:102
prepare data for inner particle dynamics
Definition: base_particle_dynamics.h:216
The base relation between a SPH body and its contact SPH bodies.
Definition: base_body_relation.h:136
A place to put all global variables.
Definition: base_particle_dynamics.h:128
Definition: base_particle_dynamics.h:77
Definition: base_particle_dynamics.h:92
std::function< ReturnType(size_t, Real)> ReduceFunctor
Definition: base_particle_dynamics.h:52
BaseMaterial * base_material_
Definition: base_body.h:87
Base of all materials.
Definition: base_material.h:52
There are the classes for neighboring particles. It saves the information for carring out pair intera...
void ParticleIterator(size_t total_real_particles, const ParticleFunctor &particle_functor, Real dt)
Definition: base_particle_dynamics.cpp:12
Definition: base_particle_dynamics.h:82
This is the base classes of SPH bodies. The real body is for that with cell linked list and the ficti...
ReturnType ReduceIterator(size_t total_real_particles, ReturnType temp, ReduceFunctor< ReturnType > &reduce_functor, ReduceOperation &reduce_operation, Real dt=0.0)
Definition: base_particle_dynamics.hpp:43
StdVec< ParticleConfiguration * > contact_configuration_
Definition: base_particle_dynamics.h:251
StdLargeVec< Neighborhood > ParticleConfiguration
Definition: neighbor_relation.h:66
prepare data for complex particle dynamics
Definition: base_particle_dynamics.h:264
void ParticleIterator_parallel(size_t total_real_particles, const ParticleFunctor &particle_functor, Real dt)
Definition: base_particle_dynamics.cpp:18
virtual void setupDynamics(Real dt=0.0)
Definition: base_particle_dynamics.h:164
Definition: base_particle_dynamics.h:113
prepare data for contact particle dynamics
Definition: base_particle_dynamics.h:240
The relation combined an inner and a contact body relation. The interaction is in a inner-boundary-co...
Definition: complex_body_relation.h:42
void ParticleIteratorSplittingSweep_parallel(SplitCellLists &split_cell_lists, const ParticleFunctor &particle_functor, Real dt)
Definition: base_particle_dynamics.cpp:61
Set up of basic data structure.
void ParticleIteratorSplittingSweep(SplitCellLists &split_cell_lists, const ParticleFunctor &particle_functor, Real dt)
Definition: base_particle_dynamics.cpp:28
BaseParticles * base_particles_
Definition: base_body.h:88
SPHBody is a base body with basic data and functions. Its derived class can be a real fluid body...
Definition: base_body.h:61
ParticleConfiguration inner_configuration_
Definition: base_body_relation.h:124
The abstract relation within a SPH body.
Definition: base_body_relation.h:117
prepare data for simple particle dynamics.
Definition: base_particle_dynamics.h:185
StdLargeVec< size_t > unsorted_id_
Definition: base_particles.h:156
Definition: base_particle_dynamics.h:171
std::function< void(size_t, Real)> ParticleFunctor
Definition: base_particle_dynamics.h:49
ReturnType ReduceIterator_parallel(size_t total_real_particles, ReturnType temp, ReduceFunctor< ReturnType > &reduce_functor, ReduceOperation &reduce_operation, Real dt=0.0)
Definition: base_particle_dynamics.hpp:54
StdLargeVec< size_t > sorted_id_
Definition: base_particles.h:157
Definition: base_particle_dynamics.h:97
Definition: solid_body_supplementary.cpp:9
particle dynamics by considering contribution from extra contact bodies
Definition: base_particle_dynamics.h:281
The base class for all particle dynamics This class contains the only two interface functions availab...
Definition: base_particle_dynamics.h:145