31 #ifndef BASE_PARTICLE_DYNAMICS_H 32 #define BASE_PARTICLE_DYNAMICS_H 34 #include "base_data_package.h" 37 #include "all_body_relations.h" 51 template <
class ReturnType>
60 template <
class ReturnType,
typename ReduceOperation>
61 ReturnType
ReduceIterator(
size_t total_real_particles, ReturnType temp,
64 template <
class ReturnType,
typename ReduceOperation>
76 template <
class ReturnType>
79 ReturnType operator()(
const ReturnType &x,
const ReturnType &y)
const {
return x + y; };
84 Real operator()(Real x, Real y)
const {
return SMAX(x, y); };
89 Real operator()(Real x, Real y)
const {
return SMIN(x, y); };
94 bool operator()(
bool x,
bool y)
const {
return x || y; };
99 bool operator()(
bool x,
bool y)
const {
return x && y; };
104 Vecd operator()(
const Vecd &x,
const Vecd &y)
const 107 for (
int i = 0; i < lower_bound.size(); ++i)
108 lower_bound[i] = SMIN(x[i], y[i]);
115 Vecd operator()(
const Vecd &x,
const Vecd &y)
const 118 for (
int i = 0; i < upper_bound.size(); ++i)
119 upper_bound[i] = SMAX(x[i], y[i]);
135 static Real physical_time_;
144 template <
class ReturnType =
void>
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;
162 void setBodyUpdated() { sph_body_->setNewlyUpdated(); };
182 template <
class BodyType =
SPHBody,
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_)),
196 BodyType *getBody() {
return body_; };
197 ParticlesType *getParticles() {
return particles_; };
198 MaterialType *getMaterial() {
return material_; };
202 ParticlesType *particles_;
203 MaterialType *material_;
204 StdLargeVec<size_t> &sorted_id_;
205 StdLargeVec<size_t> &unsorted_id_;
212 template <
class BodyType =
SPHBody,
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,
247 StdVec<ContactBodyType *> contact_bodies_;
248 StdVec<ContactParticlesType *> contact_particles_;
249 StdVec<ContactMaterialType *> contact_material_;
258 template <
class BodyType =
SPHBody,
261 class ContactBodyType =
SPHBody,
266 ContactBodyType, ContactParticlesType, ContactMaterialType, DataDelegateEmptyBase>
271 DataDelegateContact<BodyType, ParticlesType, MaterialType, ContactBodyType, ContactParticlesType,
280 template <
class ParticleDynamicsInnerType,
class ContactDataType>
294 virtual void prepareContactData() = 0;
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
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
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