SPHinXsys  alpha version
particle_dynamics_bodypart.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  * --------------------------------------------------------------------------*/
32 #ifndef PARTICLE_DYNAMICS_BODYPART_H
33 #define PARTICLE_DYNAMICS_BODYPART_H
34 
35 #include "base_particle_dynamics.h"
37 
38 #include "base_body_part.h"
39 namespace SPH
40 {
42  void PartIteratorByParticle(const IndexVector &body_part_particles, const ParticleFunctor &particle_functor, Real dt = 0.0);
44  void PartIteratorByParticle_parallel(const IndexVector &body_part_particles, const ParticleFunctor &particle_functor, Real dt = 0.0);
46  void PartIteratorByCell(const CellLists &body_part_cells, const ParticleFunctor &particle_functor, Real dt = 0.0);
48  void PartIteratorByCell_parallel(const CellLists &body_part_cells, const ParticleFunctor &particle_functor, Real dt = 0.0);
49 
57  {
58  public:
59  PartDynamicsByParticle(SPHBody &sph_body, BodyPartByParticle &body_part);
60  virtual ~PartDynamicsByParticle(){};
61 
62  protected:
63  IndexVector &body_part_particles_;
64  };
65 
71  {
72  public:
74  virtual ~PartSimpleDynamicsByParticle(){};
75 
76  virtual void exec(Real dt = 0.0) override;
77  virtual void parallel_exec(Real dt = 0.0) override;
78 
79  protected:
80  ParticleFunctor functor_update_;
81  virtual void Update(size_t index_i, Real dt = 0.0) = 0;
82  };
83 
89  {
90  public:
93 
94  virtual void exec(Real dt = 0.0) override;
95  virtual void parallel_exec(Real dt = 0.0) override;
96 
97  protected:
98  ParticleFunctor functor_interaction_;
99  virtual void Interaction(size_t index_i, Real dt = 0.0) = 0;
100  };
101 
107  {
108  public:
111 
112  virtual void exec(Real dt = 0.0) override;
113  virtual void parallel_exec(Real dt = 0.0) override;
114 
115  protected:
116  virtual void Update(size_t index_i, Real dt = 0.0) = 0;
117  ParticleFunctor functor_update_;
118  };
119 
125  {
126  public:
129 
130  virtual void exec(Real dt = 0.0) override;
131  virtual void parallel_exec(Real dt = 0.0) override;
132 
133  protected:
134  virtual void Initialization(size_t index_i, Real dt = 0.0) = 0;
135  ParticleFunctor functor_initialization_;
136  };
137 
144  {
145  public:
146  PartDynamicsByCell(SPHBody &sph_body, BodyPartByCell &body_part);
147  virtual ~PartDynamicsByCell(){};
148 
149  virtual void exec(Real dt = 0.0) override;
150  virtual void parallel_exec(Real dt = 0.0) override;
151 
152  protected:
153  CellLists &body_part_cells_;
154  ParticleFunctor functor_update_;
155  virtual void Update(size_t index_i, Real dt = 0.0) = 0;
156  };
161  template <class ReturnType, typename ReduceOperation>
162  class PartDynamicsByCellReduce : public ParticleDynamics<ReturnType>
163  {
164  public:
165  PartDynamicsByCellReduce(SPHBody &sph_body, BodyPartByCell &body_part)
166  : ParticleDynamics<ReturnType>(sph_body), body_part_cells_(body_part.body_part_cells_),
167  quantity_name_("ReducedQuantity"), initial_reference_(){};
168  virtual ~PartDynamicsByCellReduce(){};
169 
170  ReturnType InitialReference() { return initial_reference_; };
171  std::string QuantityName() { return quantity_name_; };
172 
173  virtual ReturnType exec(Real dt = 0.0) override
174  {
175  ReturnType temp = initial_reference_;
176  this->SetupReduce();
177  for (size_t i = 0; i != body_part_cells_.size(); ++i)
178  {
179  ListDataVector &list_data = body_part_cells_[i]->cell_list_data_;
180  for (size_t num = 0; num < list_data.size(); ++num)
181  {
182  temp = reduce_operation_(temp, ReduceFunction(list_data[num].first, dt));
183  }
184  }
185  return OutputResult(temp);
186  };
187 
188  virtual ReturnType parallel_exec(Real dt = 0.0) override
189  {
190  ReturnType temp = initial_reference_;
191  this->SetupReduce();
192  temp = parallel_reduce(
193  blocked_range<size_t>(0, body_part_cells_.size()),
194  temp,
195  [&](const blocked_range<size_t> &r, ReturnType temp0) -> ReturnType
196  {
197  for (size_t i = r.begin(); i != r.end(); ++i)
198  {
199  ListDataVector &list_data = body_part_cells_[i]->cell_list_data_;
200  for (size_t num = 0; num < list_data.size(); ++num)
201  {
202  temp0 = reduce_operation_(temp0, ReduceFunction(list_data[num].first, dt));
203  }
204  }
205  return temp0;
206  },
207  [this](ReturnType x, ReturnType y) -> ReturnType
208  { return reduce_operation_(x, y); });
209 
210  return OutputResult(temp);
211  };
212 
213  protected:
214  ReduceOperation reduce_operation_;
215  CellLists &body_part_cells_;
216  std::string quantity_name_;
217  ReturnType initial_reference_;
218  virtual void SetupReduce(){};
219  virtual ReturnType ReduceFunction(size_t index_i, Real dt = 0.0) = 0;
220  virtual ReturnType OutputResult(ReturnType reduced_value) { return reduced_value; };
221  };
226  template <class ReturnType, typename ReduceOperation>
228  {
229  public:
231  : ParticleDynamics<ReturnType>(sph_body),
232  body_part_particles_(body_part.body_part_particles_),
233  quantity_name_("ReducedQuantity"), initial_reference_(){};
234  virtual ~PartDynamicsByParticleReduce(){};
235 
236  ReturnType InitialReference() { return initial_reference_; };
237  std::string QuantityName() { return quantity_name_; };
238 
239  virtual ReturnType exec(Real dt = 0.0) override
240  {
241  ReturnType temp = initial_reference_;
242  this->SetupReduce();
243  for (size_t i = 0; i < body_part_particles_.size(); ++i)
244  {
245  temp = reduce_operation_(temp, ReduceFunction(body_part_particles_[i], dt));
246  }
247  return OutputResult(temp);
248  };
249  virtual ReturnType parallel_exec(Real dt = 0.0) override
250  {
251  ReturnType temp = initial_reference_;
252  this->SetupReduce();
253  temp = parallel_reduce(
254  blocked_range<size_t>(0, body_part_particles_.size()),
255  temp,
256  [&](const blocked_range<size_t> &r, ReturnType temp0) -> ReturnType
257  {
258  for (size_t n = r.begin(); n != r.end(); ++n)
259  {
260  temp0 = reduce_operation_(temp0, ReduceFunction(body_part_particles_[n], dt));
261  }
262  return temp0;
263  },
264  [this](ReturnType x, ReturnType y) -> ReturnType
265  {
266  return reduce_operation_(x, y);
267  });
268 
269  return OutputResult(temp);
270  };
271 
272  protected:
273  ReduceOperation reduce_operation_;
274  IndexVector &body_part_particles_;
275  std::string quantity_name_;
276  ReturnType initial_reference_;
277  virtual void SetupReduce(){};
278  virtual ReturnType ReduceFunction(size_t index_i, Real dt = 0.0) = 0;
279  virtual ReturnType OutputResult(ReturnType reduced_value) { return reduced_value; };
280  };
281 }
282 #endif // PARTICLE_DYNAMICS_BODYPART_H
Abstract class for particle interaction involving in a body part.
Definition: particle_dynamics_bodypart.h:88
A body part with a collection of cell lists.
Definition: base_body_part.h:100
StdVec< size_t > IndexVector
Definition: sph_data_containers.h:37
This is the implementation of the template class for 3D build.
This is the base classes of body parts.
A body part with a collection of particles.
Definition: base_body_part.h:64
StdLargeVec< CellList * > CellLists
Definition: sph_data_containers.h:46
void PartIteratorByCell_parallel(const CellLists &body_part_cells, const ParticleFunctor &particle_functor, Real dt)
Definition: particle_dynamics_bodypart.cpp:44
virtual ReturnType exec(Real dt=0.0) override
Definition: particle_dynamics_bodypart.h:239
Abstract class for reduce operation in a Eulerian constrain region.
Definition: particle_dynamics_bodypart.h:162
virtual void exec(Real dt=0.0) override
Definition: particle_dynamics_bodypart.cpp:88
Abstract class for body part simple particle dynamics.
Definition: particle_dynamics_bodypart.h:70
StdLargeVec< ListData > ListDataVector
Definition: sph_data_containers.h:44
IndexVector body_part_particles_
Definition: base_body_part.h:67
Abstract class for particle interaction involving in a body part with an extra update step...
Definition: particle_dynamics_bodypart.h:106
virtual ReturnType exec(Real dt=0.0) override
Definition: particle_dynamics_bodypart.h:173
This is for the base classes of particle dynamics, which describe the interaction between particles...
virtual void exec(Real dt=0.0) override
Definition: particle_dynamics_bodypart.cpp:69
CellLists body_part_cells_
Definition: base_body_part.h:103
Abstract class for imposing body part dynamics by particles. That is the constrained particles will b...
Definition: particle_dynamics_bodypart.h:56
Abstract class for imposing Eulerian constrain to a body. The constrained particles are in the tagged...
Definition: particle_dynamics_bodypart.h:143
virtual void exec(Real dt=0.0) override
Definition: particle_dynamics_bodypart.cpp:107
SPHBody is a base body with basic data and functions. Its derived class can be a real fluid body...
Definition: base_body.h:61
reduce operation in a Lagrangian contrained region.
Definition: particle_dynamics_bodypart.h:227
Definition: particle_dynamics_bodypart.h:124
virtual void exec(Real dt=0.0) override
Definition: particle_dynamics_bodypart.cpp:152
void PartIteratorByParticle(const IndexVector &body_part_particles, const ParticleFunctor &particle_functor, Real dt)
Definition: particle_dynamics_bodypart.cpp:12
void PartIteratorByParticle_parallel(const IndexVector &body_part_particles, const ParticleFunctor &particle_functor, Real dt)
Definition: particle_dynamics_bodypart.cpp:20
std::function< void(size_t, Real)> ParticleFunctor
Definition: base_particle_dynamics.h:49
virtual void exec(Real dt=0.0) override
Definition: particle_dynamics_bodypart.cpp:129
void PartIteratorByCell(const CellLists &body_part_cells, const ParticleFunctor &particle_functor, Real dt)
Definition: particle_dynamics_bodypart.cpp:34
Definition: solid_body_supplementary.cpp:9
The base class for all particle dynamics This class contains the only two interface functions availab...
Definition: base_particle_dynamics.h:145