ParticleHandler.h
Go to the documentation of this file.
1 // This file is part of the MercuryDPM project (https://www.mercurydpm.org).
2 // Copyright (c), The MercuryDPM Developers Team. All rights reserved.
3 // License: BSD 3-Clause License; see the LICENSE file in the root directory.
4 
5 #ifndef PARTICLE_HANDLER_H
6 #define PARTICLE_HANDLER_H
7 
8 #include "BaseHandler.h"
10 #include "Math/PSD.h"
11 
12 class SpeciesHandler;
13 
14 class BaseSpecies;
15 
17 {
18  MIN = 0, MAX = 1
19 };
20 
27 class ParticleHandler : public BaseHandler<BaseParticle>
28 {
29 public:
34 
39 
44 
48  ~ParticleHandler() override;
49 
53  void addExistingObject(BaseParticle* P) override;
54 
58  void addObject(BaseParticle* P) override;
59 
63  void addObject(int fromProcessor, BaseParticle* P);
64 
68  void addGhostObject(int fromProcessor, int toProcessor, BaseParticle* p);
69 
73  void addGhostObject(BaseParticle* P) override;
74 
78  void removeObject(unsigned int index) override;
79 
84  void removeGhostObject(unsigned int index);
85 
89  void removeLastObject();
90 
95 
100 
105 
111 
116 
122 
127 
132 
137 
142 
147 
149 
150  Mdouble getKineticEnergy() const;
151 
153 
154  Mdouble getMass() const;
155 
156  Vec3D getMassTimesPosition() const;
157 
158  Vec3D getCentreOfMass() const;
159 
160  Vec3D getMomentum() const;
161 
162  Vec3D getAngularMomentum() const;
163 
164  Mdouble getVolume() const;
165 
166  Mdouble getMeanRadius() const;
167 
172 
177 
182 
187 
192 
197 
202 
207 
211  //Mdouble getParticleAttributeLocal(std::function<Mdouble (BaseParticle*)> attribute, AttributeType type) const;
212 
222  template<typename T>
224  getParticleAttributeLocal(std::function<T(BaseParticle*)> attribute, AttributeType type) const
225  {
226  T attributeParticle;
227  T attributeFinal;
228 
229  if ((*this).getSize() == 0)
230  {
231  logger(WARN, "ParticleHandler is empty: returning 0.0");
232  attributeFinal = 0;
233  return 0.0;
234  }
235  else
236  {
237  //Initialise with the first particle found
238  attributeParticle = attribute(objects_[0]);
239  attributeFinal = attributeParticle;
240 
241  //Find the final attribute
242  for (BaseParticle* particle : (*this))
243  {
244  //Obtain the attribute
245  if (!(particle->isMPIParticle() || particle->isPeriodicGhostParticle() || particle->isFixed()))
246  {
247  attributeParticle = attribute(particle);
248  }
249 
250  //Decide what to do with the magnitude of the attribute
251  switch (type)
252  {
253  case AttributeType::MIN :
254  if (attributeParticle < attributeFinal)
255  {
256  attributeFinal = attributeParticle;
257  }
258  break;
259  case AttributeType::MAX :
260  if (attributeParticle > attributeFinal)
261  {
262  attributeFinal = attributeParticle;
263  }
264  break;
265  default :
266  logger(ERROR, "Attribute type is not recognised");
267  break;
268  }
269  }
270  return attributeFinal;
271  }
272  }
273 
283  template<typename T>
285  getParticleAttribute(std::function<T(BaseParticle*)> attribute, AttributeType type) const
286  {
287 #ifdef MERCURYDPM_USE_MPI
288  T particleAttributeLocal = getParticleAttributeLocal(attribute, type);
289  T particleAttributeGlobal;
290 
291  //Communicate local to global using the approriate rule
292  MPIContainer& communicator = MPIContainer::Instance();
293  switch(type)
294  {
295  case AttributeType::MIN :
296  communicator.allReduce(particleAttributeLocal,particleAttributeGlobal,MPI_MIN);
297  break;
298  case AttributeType::MAX :
299  communicator.allReduce(particleAttributeLocal,particleAttributeGlobal,MPI_MAX);
300  break;
301  default :
302  logger(ERROR,"Attribute type is not recognised");
303  break;
304  }
305  return particleAttributeGlobal;
306 #else
307  return getParticleAttributeLocal(attribute, type);
308 #endif
309  }
310 
315 
319  void clear() override;
320 
325  unsigned int getNumberOfFixedParticles() const;
326 
330  unsigned int getNumberOfUnfixedParticles() const;
331 
335  static BaseParticle* createObject(const std::string& type);
336 
340  BaseParticle* readAndCreateObject(std::istream& is);
341 
342 // /*!
343 // * \brief Create a new particle, based on the information from old-style restart data.
344 // */
345 // void readAndCreateOldObject(std::istream &is, const std::string& type);
346 
350  void readAndAddObject(std::istream& is) override;
351 
352  void write(std::ostream& os) const;
353 
357  void checkExtrema(BaseParticle* P);
358 
363 
367  void computeAllMasses(unsigned int indSpecies);
368 
372  void computeAllMasses();
373 
378  void addedFixedParticle();
379 
384  void removedFixedParticle();
385 
389  std::string getName() const override;
390 
394  unsigned int getNumberOfRealObjects() const;
395 
400  unsigned int getNumberOfObjects() const override;
401 
405  unsigned int getNumberOfFixedObjectsLocal() const;
406 
410  unsigned int getNumberOfFixedObjects() const;
411 
415  unsigned int getNumberOfRealObjectsLocal() const;
416 
417  void actionsBeforeTimeStep();
418 
419  void actionsAfterTimeStep();
420 
421  double getLiquidFilmVolume() const;
422 
426  void saveNumberPSDtoCSV(std::string csvFileName, std::vector<double> diameterBins = {}) const;
427 
431  PSD getPSD(std::vector<Mdouble> bins = {}, Mdouble scaleFactor = 1.0) const;
432 
433 private:
434 
436 
438 
439  Mdouble getMassLocal() const;
440 
441  Mdouble getVolumeLocal() const;
442 
444 
445  Mdouble getSumRadiusLocal() const;
446 
457 
462 
466  unsigned int NFixedParticles_;
467 };
468 
469 #endif
470 
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Eigen::Triplet< double > T
Definition: EigenUnitTest.cpp:11
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ WARN
@ ERROR
AttributeType
Definition: ParticleHandler.h:17
@ MIN
Definition: ParticleHandler.h:18
@ MAX
Definition: ParticleHandler.h:18
float * p
Definition: Tutorial_Map_using.cpp:9
Container to store the pointers to all objects that one creates in a simulation.
Definition: BaseHandler.h:30
std::vector< BaseParticle * > objects_
The actual list of Object pointers.
Definition: BaseHandler.h:283
Definition: BaseParticle.h:33
BaseSpecies is the class from which all other species are derived.
Definition: BaseSpecies.h:29
This class contains all information and functions required for communication between processors.
Definition: MpiContainer.h:109
static MPIContainer & Instance()
fetch the instance to be used for communication
Definition: MpiContainer.h:113
Contains a vector with radii and probabilities of a user defined particle size distribution (PSD)
Definition: PSD.h:47
Container to store pointers to all particles.
Definition: ParticleHandler.h:28
BaseParticle * getSmallestParticleLocal() const
Gets a pointer to the smallest BaseParticle (by interactionRadius) in this ParticleHandler of the loc...
Definition: ParticleHandler.cc:496
BaseParticle * getHighestPositionComponentParticle(int i) const
Gets a pointer to the particle with the highest coordinates in direction i in this ParticleHandler.
Definition: ParticleHandler.cc:886
~ParticleHandler() override
Destructor, it destructs the ParticleHandler and all BaseParticle it contains.
Definition: ParticleHandler.cc:86
std::enable_if< std::is_scalar< T >::value, T >::type getParticleAttribute(std::function< T(BaseParticle *)> attribute, AttributeType type) const
Computes an attribute type (min/max/..) of a particle attribute(position/velocity) in the global doma...
Definition: ParticleHandler.h:285
Mdouble getLargestInteractionRadiusLocal() const
Returns the largest interaction radius of the current domain.
Definition: ParticleHandler.cc:751
std::enable_if< std::is_scalar< T >::value, T >::type getParticleAttributeLocal(std::function< T(BaseParticle *)> attribute, AttributeType type) const
Computes an attribute type (min/max/..) of a particle attribute (position/velocity) in a local domain...
Definition: ParticleHandler.h:224
BaseParticle * getFastestParticle() const
Definition: ParticleHandler.cc:704
Mdouble getRotationalEnergy() const
Definition: ParticleHandler.cc:580
void clear() override
Empties the whole ParticleHandler by removing all BaseParticle.
Definition: ParticleHandler.cc:971
void removeGhostObject(unsigned int index)
Removes a BaseParticle from the ParticleHandler without a global check, this is only to be done for m...
Definition: ParticleHandler.cc:416
BaseParticle * getFastestParticleLocal() const
Gets a pointer to the fastest BaseParticle in this ParticleHandler.
Definition: ParticleHandler.cc:681
Vec3D getMomentum() const
Definition: ParticleHandler.cc:660
double getLiquidFilmVolume() const
Definition: ParticleHandler.cc:1386
BaseParticle * getHighestVelocityComponentParticle(int i) const
Gets a pointer to the particle with the highest velocity in direction i in this ParticleHandler.
Definition: ParticleHandler.cc:958
void readAndAddObject(std::istream &is) override
Definition: ParticleHandler.cc:1092
Mdouble getHighestPositionX() const
Function returns the highest position in the x-direction.
Definition: ParticleHandler.cc:989
static BaseParticle * createObject(const std::string &type)
Reads BaseParticle into the ParticleHandler from restart data.
Definition: ParticleHandler.cc:1019
Mdouble getKineticEnergyLocal() const
Definition: ParticleHandler.cc:538
void computeSmallestParticle()
Computes the smallest particle (by interaction radius) and sets it in smallestParticle_.
Definition: ParticleHandler.cc:445
Mdouble getSmallestInteractionRadiusLocal() const
Returns the smallest interaction radius of the current domain.
Definition: ParticleHandler.cc:715
std::string getName() const override
Returns the name of the handler, namely the string "ParticleHandler".
Definition: ParticleHandler.cc:1245
Vec3D getMassTimesPositionLocal() const
Definition: ParticleHandler.cc:621
Vec3D getMassTimesPosition() const
Definition: ParticleHandler.cc:630
Mdouble getMassLocal() const
Definition: ParticleHandler.cc:596
Vec3D getCentreOfMass() const
Definition: ParticleHandler.cc:648
BaseParticle * getLowestVelocityComponentParticleLocal(int i) const
Gets a pointer to the particle with the lowest velocity in direction i in this ParticleHandler.
Definition: ParticleHandler.cc:899
Vec3D getAngularMomentum() const
Definition: ParticleHandler.cc:669
unsigned int getNumberOfRealObjectsLocal() const
Returns the number of real objects on a local domain. MPI particles and periodic particles are neglec...
Definition: ParticleHandler.cc:1281
BaseParticle * getLowestPositionComponentParticleLocal(int i) const
Gets a pointer to the particle with the lowest coordinates in direction i in this ParticleHandler.
Definition: ParticleHandler.cc:823
void actionsAfterTimeStep()
Definition: ParticleHandler.cc:1378
BaseParticle * getLowestVelocityComponentParticle(int i) const
Gets a pointer to the particle with the lowest velocity in direction i in this ParticleHandler.
Definition: ParticleHandler.cc:922
void addExistingObject(BaseParticle *P) override
Adds a BaseParticle to the ParticleHandler.
Definition: ParticleHandler.cc:110
void addObject(BaseParticle *P) override
Adds a BaseParticle to the ParticleHandler.
Definition: ParticleHandler.cc:150
ParticleHandler & operator=(const ParticleHandler &rhs)
Assignment operator.
Definition: ParticleHandler.cc:64
Mdouble getLargestInteractionRadius() const
Returns the largest interaction radius.
Definition: ParticleHandler.cc:766
Mdouble getSmallestInteractionRadius() const
Returns the smallest interaction radius.
Definition: ParticleHandler.cc:730
unsigned int getNumberOfFixedObjects() const
Computes the number of fixed particles in the whole simulation.
Definition: ParticleHandler.cc:1356
void write(std::ostream &os) const
BaseParticle * readAndCreateObject(std::istream &is)
Create a new particle, based on the information provided in a restart file.
Definition: ParticleHandler.cc:1063
unsigned int getNumberOfRealObjects() const
Returns the number of real objects (on all processors)
Definition: ParticleHandler.cc:1302
void actionsBeforeTimeStep()
Definition: ParticleHandler.cc:1370
Mdouble getVolume() const
Definition: ParticleHandler.cc:1261
BaseParticle * smallestParticle_
A pointer to the smallest BaseParticle (by interactionRadius) in this ParticleHandler.
Definition: ParticleHandler.h:461
void checkExtrema(BaseParticle *P)
Checks if the extrema of this ParticleHandler needs updating.
Definition: ParticleHandler.cc:1171
void removeObject(unsigned int index) override
Removes a BaseParticle from the ParticleHandler.
Definition: ParticleHandler.cc:388
unsigned int NFixedParticles_
Number of fixed particles.
Definition: ParticleHandler.h:466
BaseParticle * getLargestParticleLocal() const
Gets a pointer to the largest BaseParticle (by interactionRadius) in the ParticleHandler of the local...
Definition: ParticleHandler.cc:520
unsigned int getNumberOfObjects() const override
Returns the number of objects in the container. In parallel code this practice is forbidden to avoid ...
Definition: ParticleHandler.cc:1323
BaseParticle * largestParticle_
A pointer to the largest BaseParticle (by interactionRadius) in this ParticleHandler.
Definition: ParticleHandler.h:456
unsigned int getNumberOfFixedObjectsLocal() const
Computes the number of Fixed particles on a local domain.
Definition: ParticleHandler.cc:1340
BaseParticle * getHighestPositionComponentParticleLocal(int i) const
Gets a pointer to the particle with the highest coordinates in direction i in this ParticleHandler.
Definition: ParticleHandler.cc:862
Mdouble getVolumeLocal() const
Definition: ParticleHandler.cc:1250
void checkExtremaOnDelete(BaseParticle *P)
Checks if the extrema of this ParticleHandler needs updating when a particle is deleted.
Definition: ParticleHandler.cc:1197
Mdouble getKineticEnergy() const
Definition: ParticleHandler.cc:551
void addGhostObject(int fromProcessor, int toProcessor, BaseParticle *p)
Adds a ghost particle located at fromProcessor to toProcessor.
Definition: ParticleHandler.cc:275
Mdouble getSumRadiusLocal() const
Definition: ParticleHandler.cc:785
BaseParticle * getLowestPositionComponentParticle(int i) const
Gets a pointer to the particle with the lowest coordinates in direction i in this ParticleHandler.
Definition: ParticleHandler.cc:849
void removeLastObject()
Removes the last BaseParticle from the ParticleHandler.
Definition: ParticleHandler.cc:435
BaseParticle * getHighestVelocityComponentParticleLocal(int i) const
Gets a pointer to the particle with the highest velocity in direction i in this ParticleHandler.
Definition: ParticleHandler.cc:935
Mdouble getRotationalEnergyLocal() const
Definition: ParticleHandler.cc:567
unsigned int getNumberOfFixedParticles() const
Gets the number of particles that are fixed.
Definition: ParticleHandler.cc:1003
void computeAllMasses()
Computes the mass for all BaseParticle in this ParticleHandler.
Definition: ParticleHandler.cc:1224
BaseParticle * getLargestParticle() const
Returns the pointer of the largest particle in the particle handler. When mercury is running in paral...
Definition: ParticleHandler.cc:528
BaseParticle * getSmallestParticle() const
Gets a pointer to the smallest BaseParticle (by interactionRadius) in this ParticleHandler of the loc...
Definition: ParticleHandler.cc:505
PSD getPSD(std::vector< Mdouble > bins={}, Mdouble scaleFactor=1.0) const
Gets the PSD from the particles currently in the handler.
Definition: ParticleHandler.cc:1427
ParticleHandler()
Default constructor, it creates an empty ParticleHandler.
Definition: ParticleHandler.cc:29
unsigned int getNumberOfUnfixedParticles() const
Gets the number of particles that are not fixed.
Definition: ParticleHandler.cc:1011
void addedFixedParticle()
Increment of the number of fixed particles.
Definition: ParticleHandler.cc:1232
void computeLargestParticle()
Computes the largest particle (by interaction radius) and sets it in largestParticle_.
Definition: ParticleHandler.cc:469
void removedFixedParticle()
Decrement of the number of fixed particles.
Definition: ParticleHandler.cc:1237
Mdouble getMeanRadius() const
Definition: ParticleHandler.cc:798
Mdouble getMass() const
Definition: ParticleHandler.cc:605
void saveNumberPSDtoCSV(std::string csvFileName, std::vector< double > diameterBins={}) const
Writes the PSD of the particles currently in the handler to a CSV file, with type PROBABILITYDENSITY_...
Definition: ParticleHandler.cc:1400
Container to store all ParticleSpecies.
Definition: SpeciesHandler.h:15
Definition: Kernel/Math/Vector.h:30
squared absolute value
Definition: GlobalFunctions.h:87
double P
Uniform pressure.
Definition: TwenteMeshGluing.cpp:77
bins
Definition: UniformPSDSelfTest.py:19
default
Definition: calibrate.py:45
type
Definition: compute_granudrum_aor.py:141
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286