SPHinXsys  alpha version
relax_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  * --------------------------------------------------------------------------*/
30 #ifndef RELAX_DYNAMICS_H
31 #define RELAX_DYNAMICS_H
32 
33 #include "all_particle_dynamics.h"
34 #include "base_kernel.h"
35 #include "cell_linked_list.h"
36 #include "all_body_relations.h"
37 #include "general_dynamics.h"
38 
39 namespace SPH
40 {
41  class GeometryShape;
42  class LevelSetShape;
43 
44  namespace relax_dynamics
45  {
46  typedef DataDelegateSimple<SPHBody, BaseParticles> RelaxDataDelegateSimple;
47 
48  typedef DataDelegateInner<SPHBody, BaseParticles> RelaxDataDelegateInner;
49 
50  typedef DataDelegateComplex<SPHBody, BaseParticles, BaseMaterial, SPHBody, BaseParticles> RelaxDataDelegateComplex;
51 
57  class GetTimeStepSizeSquare : public ParticleDynamicsReduce<Real, ReduceMax>,
59  {
60  public:
61  explicit GetTimeStepSizeSquare(SPHBody &sph_body);
62  virtual ~GetTimeStepSizeSquare(){};
63 
64  protected:
65  StdLargeVec<Vecd> &acc_;
66  Real h_ref_;
67  Real ReduceFunction(size_t index_i, Real dt = 0.0) override;
68  Real OutputResult(Real reduced_value) override;
69  };
70 
78  {
79  public:
80  explicit RelaxationAccelerationInner(BaseBodyRelationInner &inner_relation);
81  virtual ~RelaxationAccelerationInner(){};
82 
83  protected:
84  StdLargeVec<Real> &Vol_;
85  StdLargeVec<Vecd> &acc_;
86  StdLargeVec<Vecd> &pos_;
87  virtual void Interaction(size_t index_i, Real dt = 0.0) override;
88  };
89 
95  {
96  public:
98  BaseBodyRelationInner &inner_relation);
100 
101  protected:
102  LevelSetShape *level_set_shape_;
103  virtual void Interaction(size_t index_i, Real dt = 0.0) override;
104  };
105 
112  {
113  public:
114  explicit UpdateParticlePosition(SPHBody &sph_body);
115  virtual ~UpdateParticlePosition(){};
116 
117  protected:
118  StdLargeVec<Vecd> &pos_, &acc_;
119  virtual void Update(size_t index_i, Real dt = 0.0) override;
120  };
121 
128  {
129  public:
132 
133  protected:
134  StdLargeVec<Real> &h_ratio_, &Vol_;
135  StdLargeVec<Vecd> &pos_;
136  Shape &body_shape_;
137  Kernel &kernel_;
138  ParticleSpacingByBodyShape *particle_spacing_by_body_shape_;
139 
140  virtual void Update(size_t index_i, Real dt = 0.0) override;
141  };
142 
151  {
152  public:
153  explicit RelaxationAccelerationComplex(ComplexBodyRelation &body_complex_relation);
154  virtual ~RelaxationAccelerationComplex(){};
155 
156  protected:
157  StdLargeVec<Real> &Vol_;
158  StdLargeVec<Vecd> &acc_, &pos_;
159  StdVec<StdLargeVec<Real> *> contact_Vol_;
160  virtual void Interaction(size_t index_i, Real dt = 0.0) override;
161  };
162 
171  {
172  public:
173  ShapeSurfaceBounding(SPHBody &sph_body, NearShapeSurface &body_part);
174  virtual ~ShapeSurfaceBounding(){};
175 
176  protected:
177  StdLargeVec<Vecd> &pos_;
178  LevelSetShape *level_set_shape_;
179  Real constrained_distance_;
180  virtual void Update(size_t index_i, Real dt = 0.0) override;
181  };
182 
191  {
192  public:
193  ConstraintSurfaceParticles(SPHBody &sph_body, BodySurface &body_part);
194  virtual ~ConstraintSurfaceParticles(){};
195 
196  protected:
197  Real constrained_distance_;
198  StdLargeVec<Vecd> &pos_;
199  LevelSetShape *level_set_shape_;
200  virtual void Update(size_t index_i, Real dt = 0.0) override;
201  };
202 
208  {
209  protected:
210  RealBody *real_body_;
211  BaseBodyRelationInner &inner_relation_;
212  NearShapeSurface near_shape_surface_;
213 
214  public:
215  explicit RelaxationStepInner(BaseBodyRelationInner &inner_relation,
216  bool level_set_correction = false);
217  virtual ~RelaxationStepInner(){};
218 
219  UniquePtr<RelaxationAccelerationInner> relaxation_acceleration_inner_;
220  GetTimeStepSizeSquare get_time_step_square_;
221  UpdateParticlePosition update_particle_position_;
222  ShapeSurfaceBounding surface_bounding_;
223 
224  virtual void exec(Real dt = 0.0) override;
225  virtual void parallel_exec(Real dt = 0.0) override;
226  };
227 
236  {
237  public:
239  ComplexBodyRelation &body_complex_relation, const std::string &shape_name);
241 
242  protected:
243  LevelSetShape *level_set_shape_;
244  virtual void Interaction(size_t index_i, Real dt = 0.0) override;
245  };
246 
252  {
253  protected:
254  RealBody *real_body_;
255  ComplexBodyRelation &complex_relation_;
256  NearShapeSurface near_shape_surface_;
257 
258  public:
259  explicit RelaxationStepComplex(ComplexBodyRelation &body_complex_relation,
260  const std::string &shape_name, bool level_set_correction = false);
261  virtual ~RelaxationStepComplex(){};
262 
263  UniquePtr<RelaxationAccelerationComplex> relaxation_acceleration_complex_;
264  GetTimeStepSizeSquare get_time_step_square_;
265  UpdateParticlePosition update_particle_position_;
266  ShapeSurfaceBounding surface_bounding_;
267 
268  virtual void exec(Real dt = 0.0) override;
269  virtual void parallel_exec(Real dt = 0.0) override;
270  };
271 
281  {
282  public:
283  ShellMidSurfaceBounding(SPHBody &body, NearShapeSurface &body_part, BaseBodyRelationInner &inner_relation,
284  Real thickness, Real level_set_refinement_ratio);
285  virtual ~ShellMidSurfaceBounding(){};
286 
287  protected:
288  StdLargeVec<Vecd> &pos_;
289  Real constrained_distance_;
290  LevelSetShape *level_set_shape_;
291  Real particle_spacing_ref_, thickness_, level_set_refinement_ratio_;
292  virtual void Update(size_t index_i, Real dt = 0.0) override;
293  };
294 
300  {
301  const Real convergence_criterion_;
302  const Real consistency_criterion_;
303 
304  void predictNormalDirection();
305  void correctNormalDirection();
306 
307  public:
308  explicit ShellNormalDirectionPrediction(BaseBodyRelationInner &inner_relation,
309  Real thickness, Real consistency_criterion = cos(Pi / 20.0));
310  virtual ~ShellNormalDirectionPrediction(){};
311 
312  virtual void exec(Real dt = 0.0) override;
313  virtual void parallel_exec(Real dt = 0.0) override { exec(); };
314 
315  protected:
317  {
318  Real thickness_;
319  LevelSetShape *level_set_shape_;
320  StdLargeVec<Vecd> &pos_, &n_, n_temp_;
321 
322  public:
323  NormalPrediction(SPHBody &sph_body, Real thickness);
324  virtual ~NormalPrediction(){};
325  void update(size_t index_i, Real dt = 0.0);
326  };
327 
328  class PredictionConvergenceCheck : public ParticleDynamicsReduce<bool, ReduceAND>,
330  {
331  public:
332  PredictionConvergenceCheck(SPHBody &sph_body, Real convergence_criterion);
333  virtual ~PredictionConvergenceCheck(){};
334 
335  protected:
336  const Real convergence_criterion_;
337  StdLargeVec<Vecd> &n_, &n_temp_;
338  bool ReduceFunction(size_t index_i, Real dt = 0.0) override;
339  };
340 
342  {
343  public:
344  explicit ConsistencyCorrection(BaseBodyRelationInner &inner_relation, Real consistency_criterion);
345  virtual ~ConsistencyCorrection(){};
346 
348  virtual void parallel_exec(Real dt = 0.0) override { exec(); };
349 
350  protected:
351  const Real consistency_criterion_;
352  StdLargeVec<int> updated_indicator_;
354  virtual void Interaction(size_t index_i, Real dt = 0.0) override;
355  };
356 
357  class ConsistencyUpdatedCheck : public ParticleDynamicsReduce<bool, ReduceAND>,
359  {
360  public:
361  explicit ConsistencyUpdatedCheck(SPHBody &sph_body);
362  virtual ~ConsistencyUpdatedCheck(){};
363 
364  protected:
365  StdLargeVec<int> &updated_indicator_;
366  bool ReduceFunction(size_t index_i, Real dt = 0.0) override;
367  };
368 
369  class SmoothingNormal : public ParticleSmoothing<Vecd>
370  {
371  public:
372  explicit SmoothingNormal(BaseBodyRelationInner &inner_relation);
373  virtual ~SmoothingNormal(){};
374 
375  protected:
376  virtual void Update(size_t index_i, Real dt = 0.0) override;
377  };
378 
379  SimpleDynamics<NormalPrediction> normal_prediction_;
380  PredictionConvergenceCheck normal_prediction_convergence_check_;
381  ConsistencyCorrection consistency_correction_;
382  ConsistencyUpdatedCheck consistency_updated_check_;
383  SmoothingNormal smoothing_normal_;
384  };
385 
391  {
392  public:
393  explicit ShellRelaxationStepInner(BaseBodyRelationInner &inner_relation, Real thickness,
394  Real level_set_refinement_ratio, bool level_set_correction = false);
395  virtual ~ShellRelaxationStepInner(){};
396 
397  UpdateParticlePosition update_shell_particle_position_;
398  ShellMidSurfaceBounding mid_surface_bounding_;
399 
400  virtual void exec(Real dt = 0.0) override;
401  virtual void parallel_exec(Real dt = 0.0) override;
402  };
403  }
404 }
405 #endif // RELAX_DYNAMICS_H
computing smoothed variable field by averaging with neighbors
Definition: general_dynamics.h:95
compute relaxation acceleration while consider the present of contact bodies with considering contact...
Definition: relax_dynamics.h:235
prepare data for inner particle dynamics
Definition: base_particle_dynamics.h:216
A body part with the collection of particles at surface of a body.
Definition: base_body_part.h:138
carry out particle relaxation step of particles within the shell body
Definition: relax_dynamics.h:390
Base class for all volumetric geometries Note that checkContain and findClosest point are basic funct...
Definition: base_geometry.h:64
Simple particle dynamics without considering particle interaction.
Definition: particle_dynamics_algorithms.h:48
virtual void exec(Real dt=0.0) override
Definition: relax_dynamics.cpp:238
This is the particle dynamics aplliable for all type bodies.
virtual void Interaction(size_t index_i, Real dt=0.0) override
Definition: relax_dynamics.cpp:107
Here gives the classes for managing cell linked lists. This is the basic class for building the parti...
relaxation dynamics for particle initialization computing the square of time step size ...
Definition: relax_dynamics.h:57
Abstract class for body part simple particle dynamics.
Definition: particle_dynamics_bodypart.h:70
compute relaxation acceleration while consider the present of contact bodies with considering contact...
Definition: relax_dynamics.h:149
Definition: particle_dynamics_algorithms.h:247
Abstract base class of a general SPH kernel function which is a smoothed Dirac delta function...
Definition: base_kernel.h:63
Base abstract class for reduce.
Definition: particle_dynamics_algorithms.h:69
prepare data for complex particle dynamics
Definition: base_particle_dynamics.h:264
simple algorithm for physics relaxation without considering contact interaction. this is usually used...
Definition: relax_dynamics.h:77
update the particle position for a time step
Definition: relax_dynamics.h:110
we constrain particles to a level function representing the interafce.
Definition: relax_dynamics.h:94
virtual void parallel_exec(Real dt=0.0) override
Definition: relax_dynamics.h:348
virtual void exec(Real dt=0.0) override
Definition: relax_dynamics.cpp:287
A shape using level set to define geometry.
Definition: level_set_shape.h:46
A body part with the cell lists near the surface of a prescribed shape.
Definition: base_body_part.h:187
Adaptive resolutions within a SPH body according to the distance to the body surface.
Definition: adaptation.h:128
prodict the normal direction of shell particles.
Definition: relax_dynamics.h:299
virtual void exec(Real dt=0.0) override
Definition: relax_dynamics.cpp:440
constrain surface particles by map contrained particles to geometry face and r = r + phi * norm (vect...
Definition: relax_dynamics.h:169
constrain particles by constraining particles to mid-surface. Note that level_set_refinement_ratio sh...
Definition: relax_dynamics.h:279
virtual void exec(Real dt=0.0) override
Definition: relax_dynamics.cpp:182
update the particle smoothing length ratio
Definition: relax_dynamics.h:126
The relation combined an inner and a contact body relation. The interaction is in a inner-boundary-co...
Definition: complex_body_relation.h:42
carry out particle relaxation step of particles within multi bodies
Definition: relax_dynamics.h:251
This is the class for particle interaction with other particles.
Definition: particle_dynamics_algorithms.h:115
Abstract class for imposing Eulerian constrain to a body. The constrained particles are in the tagged...
Definition: particle_dynamics_bodypart.h:143
SPHBody is a base body with basic data and functions. Its derived class can be a real fluid body...
Definition: base_body.h:61
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
This is the base classes of kernel functions. Implementation will be implemented in derived classes...
StdLargeVec< Vecd > & n_
Definition: relax_dynamics.h:353
Derived body with inner particle configuration or inner interactions. After construction, the particle and material must be specified.
Definition: base_body.h:182
carry out particle relaxation step of particles within the body
Definition: relax_dynamics.h:207
Definition: solid_body_supplementary.cpp:9
virtual void exec(Real dt=0.0) override
Definition: particle_dynamics_algorithms.cpp:31
constrain surface particles by map contrained particles to geometry face and r = r + phi * norm (vect...
Definition: relax_dynamics.h:189
The base class for all particle dynamics This class contains the only two interface functions availab...
Definition: base_particle_dynamics.h:145