ParticleHandler Class Reference

Container to store pointers to all particles. More...

#include <ParticleHandler.h>

+ Inheritance diagram for ParticleHandler:

Public Member Functions

 ParticleHandler ()
 Default constructor, it creates an empty ParticleHandler. More...
 
 ParticleHandler (const ParticleHandler &PH)
 Constructor that copies all BaseParticle it contains and then sets the smallest and largest particle. More...
 
ParticleHandleroperator= (const ParticleHandler &rhs)
 Assignment operator. More...
 
 ~ParticleHandler () override
 Destructor, it destructs the ParticleHandler and all BaseParticle it contains. More...
 
void addExistingObject (BaseParticle *P) override
 Adds a BaseParticle to the ParticleHandler. More...
 
void addObject (BaseParticle *P) override
 Adds a BaseParticle to the ParticleHandler. More...
 
void addObject (int fromProcessor, BaseParticle *P)
 Adds a BaseParticle located at processor fromProcessor to toProcessor. More...
 
void addGhostObject (int fromProcessor, int toProcessor, BaseParticle *p)
 Adds a ghost particle located at fromProcessor to toProcessor. More...
 
void addGhostObject (BaseParticle *P) override
 Adds a BaseParticle to the ParticleHandler. More...
 
void removeObject (unsigned int index) override
 Removes a BaseParticle from the ParticleHandler. More...
 
void removeGhostObject (unsigned int index)
 Removes a BaseParticle from the ParticleHandler without a global check, this is only to be done for mpi routines. More...
 
void removeLastObject ()
 Removes the last BaseParticle from the ParticleHandler. More...
 
void computeSmallestParticle ()
 Computes the smallest particle (by interaction radius) and sets it in smallestParticle_. More...
 
void computeLargestParticle ()
 Computes the largest particle (by interaction radius) and sets it in largestParticle_. More...
 
BaseParticlegetSmallestParticleLocal () const
 Gets a pointer to the smallest BaseParticle (by interactionRadius) in this ParticleHandler of the local domain. More...
 
BaseParticlegetSmallestParticle () const
 Gets a pointer to the smallest BaseParticle (by interactionRadius) in this ParticleHandler of the local domain. When mercury is running in parallel this function throws an error since a pointer to another domain is useless. More...
 
BaseParticlegetLargestParticleLocal () const
 Gets a pointer to the largest BaseParticle (by interactionRadius) in the ParticleHandler of the local domain. More...
 
BaseParticlegetLargestParticle () const
 Returns the pointer of the largest particle in the particle handler. When mercury is running in parallel this function throws an error since a pointer to another domain is useless. More...
 
Mdouble getSmallestInteractionRadiusLocal () const
 Returns the smallest interaction radius of the current domain. More...
 
Mdouble getSmallestInteractionRadius () const
 Returns the smallest interaction radius. More...
 
Mdouble getLargestInteractionRadiusLocal () const
 Returns the largest interaction radius of the current domain. More...
 
Mdouble getLargestInteractionRadius () const
 Returns the largest interaction radius. More...
 
BaseParticlegetFastestParticleLocal () const
 Gets a pointer to the fastest BaseParticle in this ParticleHandler. More...
 
BaseParticlegetFastestParticle () const
 
Mdouble getKineticEnergy () const
 
Mdouble getRotationalEnergy () const
 
Mdouble getMass () const
 
Vec3D getMassTimesPosition () const
 
Vec3D getCentreOfMass () const
 
Vec3D getMomentum () const
 
Vec3D getAngularMomentum () const
 
Mdouble getVolume () const
 
Mdouble getMeanRadius () const
 
BaseParticlegetLowestPositionComponentParticleLocal (int i) const
 Gets a pointer to the particle with the lowest coordinates in direction i in this ParticleHandler. More...
 
BaseParticlegetLowestPositionComponentParticle (int i) const
 Gets a pointer to the particle with the lowest coordinates in direction i in this ParticleHandler. More...
 
BaseParticlegetHighestPositionComponentParticleLocal (int i) const
 Gets a pointer to the particle with the highest coordinates in direction i in this ParticleHandler. More...
 
BaseParticlegetHighestPositionComponentParticle (int i) const
 Gets a pointer to the particle with the highest coordinates in direction i in this ParticleHandler. More...
 
BaseParticlegetLowestVelocityComponentParticleLocal (int i) const
 Gets a pointer to the particle with the lowest velocity in direction i in this ParticleHandler. More...
 
BaseParticlegetLowestVelocityComponentParticle (int i) const
 Gets a pointer to the particle with the lowest velocity in direction i in this ParticleHandler. More...
 
BaseParticlegetHighestVelocityComponentParticleLocal (int i) const
 Gets a pointer to the particle with the highest velocity in direction i in this ParticleHandler. More...
 
BaseParticlegetHighestVelocityComponentParticle (int i) const
 Gets a pointer to the particle with the highest velocity in direction i in this ParticleHandler. More...
 
template<typename T >
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. More...
 
template<typename T >
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 domain. More...
 
Mdouble getHighestPositionX () const
 Function returns the highest position in the x-direction. More...
 
void clear () override
 Empties the whole ParticleHandler by removing all BaseParticle. More...
 
unsigned int getNumberOfFixedParticles () const
 Gets the number of particles that are fixed. More...
 
unsigned int getNumberOfUnfixedParticles () const
 Gets the number of particles that are not fixed. More...
 
BaseParticlereadAndCreateObject (std::istream &is)
 Create a new particle, based on the information provided in a restart file. More...
 
void readAndAddObject (std::istream &is) override
 
void write (std::ostream &os) const
 
void checkExtrema (BaseParticle *P)
 Checks if the extrema of this ParticleHandler needs updating. More...
 
void checkExtremaOnDelete (BaseParticle *P)
 Checks if the extrema of this ParticleHandler needs updating when a particle is deleted. More...
 
void computeAllMasses (unsigned int indSpecies)
 Computes the mass for all BaseParticle of the given species in this ParticleHandler. More...
 
void computeAllMasses ()
 Computes the mass for all BaseParticle in this ParticleHandler. More...
 
void addedFixedParticle ()
 Increment of the number of fixed particles. More...
 
void removedFixedParticle ()
 Decrement of the number of fixed particles. More...
 
std::string getName () const override
 Returns the name of the handler, namely the string "ParticleHandler". More...
 
unsigned int getNumberOfRealObjects () const
 Returns the number of real objects (on all processors) More...
 
unsigned int getNumberOfObjects () const override
 Returns the number of objects in the container. In parallel code this practice is forbidden to avoid confusion with real and fake particles. If the size of the container is wished, call size() More...
 
unsigned int getNumberOfFixedObjectsLocal () const
 Computes the number of Fixed particles on a local domain. More...
 
unsigned int getNumberOfFixedObjects () const
 Computes the number of fixed particles in the whole simulation. More...
 
unsigned int getNumberOfRealObjectsLocal () const
 Returns the number of real objects on a local domain. MPI particles and periodic particles are neglected. More...
 
void actionsBeforeTimeStep ()
 
void actionsAfterTimeStep ()
 
double getLiquidFilmVolume () const
 
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_NUMBER_DISTRIBUTION. More...
 
PSD getPSD (std::vector< Mdouble > bins={}, Mdouble scaleFactor=1.0) const
 Gets the PSD from the particles currently in the handler. More...
 
- Public Member Functions inherited from BaseHandler< BaseParticle >
 BaseHandler ()
 Default BaseHandler constructor, it creates an empty BaseHandler and assigns DPMBase_ to a null pointer. More...
 
 BaseHandler (const BaseHandler< BaseParticle > &BH)
 Constructor that copies the objects of the given handler into itself and sets other variables to 0/nullptr. More...
 
virtual ~BaseHandler ()
 Destructor, it destructs the BaseHandler and all Object it contains. More...
 
void copyContentsFromOtherHandler (const BaseHandler< BaseParticle > &BH)
 Function that copies the contents (vector of pointers, maxObject_, nextId_, DPMBase_) from one handler (container) to the other. More...
 
std::enable_if<!std::is_pointer< U >::value, U * >::type copyAndAddObject (const U &object)
 Creates a copy of a Object and adds it to the BaseHandler. More...
 
std::enable_if< std::is_pointer< U >::value, U >::type copyAndAddObject (const U object)
 Creates a copy of a Object and adds it to the BaseHandler. More...
 
std::enable_if<!std::is_pointer< U >::value, U * >::type copyAndAddGhostObject (const U &object)
 Creates a copy of a Object and adds it to the BaseHandler. This is one locally for inserting mpi particles, they avoid the global check if the particle can actually be inserted, because the mpi domain already knows that is the case. More...
 
std::enable_if< std::is_pointer< U >::value, U >::type copyAndAddGhostObject (const U object)
 Creates a copy of a Object and adds it to the BaseHandler. This is one locally for inserting mpi particles, they avoid the global check if the particle can actually be inserted, because the mpi domain already knows that is the case. More...
 
virtual void addExistingObject (BaseParticle *O)
 Adds an existing object to the BaseHandler without changing the id of the object. More...
 
virtual void addObject (BaseParticle *object)
 Adds a new Object to the BaseHandler. More...
 
virtual void addGhostObject (BaseParticle *O)
 Adds a new Object to the BaseHandler. called by the to avoid increasing the id. More...
 
void removeIf (const std::function< bool(BaseParticle *)> cond)
 
virtual void removeObjects (std::vector< unsigned int > indices)
 
void removeLastObject ()
 Removes the last Object from the BaseHandler. More...
 
void read (std::istream &is)
 Reads all objects from restart data. More...
 
BaseParticlegetObjectById (const unsigned int id)
 Gets a pointer to the Object at the specified index in the BaseHandler. More...
 
std::vector< BaseParticle * > getObjectsById (const unsigned int id)
 Gets a vector of pointers to the objects with the specific id. More...
 
BaseParticlegetObject (const unsigned int id)
 Gets a pointer to the Object at the specified index in the BaseHandler.
More...
 
const BaseParticlegetObject (const unsigned int id) const
 Gets a constant pointer to the Object at the specified index in the BaseHandler. More...
 
BaseParticlegetLastObject ()
 Gets a pointer to the last Object in this BaseHandler. More...
 
const BaseParticlegetLastObject () const
 Gets a constant pointer to the last Object in this BaseHandler. More...
 
unsigned int getSize () const
 Gets the size of the particleHandler (including mpi and periodic particles) More...
 
unsigned int getStorageCapacity () const
 Gets the storage capacity of this BaseHandler. More...
 
void setStorageCapacity (const unsigned int N)
 Sets the storage capacity of this BaseHandler. More...
 
void resize (const unsigned int N, const BaseParticle &obj)
 Resizes the container to contain N elements. More...
 
const std::vector< BaseParticle * >::const_iterator begin () const
 Gets the begin of the const_iterator over all Object in this BaseHandler. More...
 
const std::vector< BaseParticle * >::iterator begin ()
 Gets the begin of the iterator over all BaseBoundary in this BaseHandler. More...
 
const std::vector< BaseParticle * >::const_iterator end () const
 Gets the end of the const_iterator over all BaseBoundary in this BaseHandler. More...
 
const std::vector< BaseParticle * >::iterator end ()
 Gets the end of the iterator over all BaseBoundary in this BaseHandler. More...
 
void setDPMBase (DPMBase *DPMBase)
 Sets the problem that is solved using this handler. More...
 
void setId (BaseParticle *object, unsigned int id)
 
void increaseId ()
 
unsigned int getNextId ()
 
void setNextId (unsigned int id)
 
DPMBasegetDPMBase ()
 Gets the problem that is solved using this handler. More...
 
DPMBasegetDPMBase () const
 Gets the problem that is solved using this handler and does not change the class. More...
 
virtual void writeVTK () const
 now empty function for writing VTK files. More...
 
unsigned getNextGroupId ()
 Should be called each time you assign a groupId. Returns the value of nextGroupId_ and increases nextGroupId_ by one. More...
 

Static Public Member Functions

static BaseParticlecreateObject (const std::string &type)
 Reads BaseParticle into the ParticleHandler from restart data. More...
 

Private Member Functions

Mdouble getKineticEnergyLocal () const
 
Mdouble getRotationalEnergyLocal () const
 
Mdouble getMassLocal () const
 
Mdouble getVolumeLocal () const
 
Vec3D getMassTimesPositionLocal () const
 
Mdouble getSumRadiusLocal () const
 

Private Attributes

BaseParticlelargestParticle_
 A pointer to the largest BaseParticle (by interactionRadius) in this ParticleHandler. More...
 
BaseParticlesmallestParticle_
 A pointer to the smallest BaseParticle (by interactionRadius) in this ParticleHandler. More...
 
unsigned int NFixedParticles_
 Number of fixed particles. More...
 

Additional Inherited Members

- Protected Attributes inherited from BaseHandler< BaseParticle >
std::vector< BaseParticle * > objects_
 The actual list of Object pointers. More...
 

Detailed Description

Container to store pointers to all particles.

The ParticleHandler is a container that stores pointers to all particles belonging to a simulation. It is implemented as a vector of pointers of type BaseParticle*.

Constructor & Destructor Documentation

◆ ParticleHandler() [1/2]

ParticleHandler::ParticleHandler ( )

Default constructor, it creates an empty ParticleHandler.

Constructor of the ParticleHandler class. It creates and empty ParticleHandler.

30 {
31  largestParticle_ = nullptr;
32  smallestParticle_ = nullptr;
33  logger(DEBUG, "ParticleHandler::ParticleHandler() finished");
34 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ DEBUG
BaseParticle * smallestParticle_
A pointer to the smallest BaseParticle (by interactionRadius) in this ParticleHandler.
Definition: ParticleHandler.h:461
BaseParticle * largestParticle_
A pointer to the largest BaseParticle (by interactionRadius) in this ParticleHandler.
Definition: ParticleHandler.h:456

References DEBUG, largestParticle_, logger, and smallestParticle_.

◆ ParticleHandler() [2/2]

ParticleHandler::ParticleHandler ( const ParticleHandler PH)

Constructor that copies all BaseParticle it contains and then sets the smallest and largest particle.

Parameters
[in]PHThe ParticleHandler that has to be copied.

This is not a copy constructor! It copies the DPMBase and all BaseParticle, and sets the other variables to 0. After that, it computes the smallest and largest particle in this handler.

44 {
45  clear();
46  setDPMBase(PH.getDPMBase());
47  largestParticle_ = nullptr;
48  smallestParticle_ = nullptr;
50  if (!objects_.empty())
51  {
54  }
55  logger(DEBUG, "ParticleHandler::ParticleHandler(const ParticleHandler &PH) finished");
56 }
void setDPMBase(DPMBase *DPMBase)
Sets the problem that is solved using this handler.
Definition: BaseHandler.h:726
std::vector< BaseParticle * > objects_
The actual list of Object pointers.
Definition: BaseHandler.h:283
DPMBase * getDPMBase()
Gets the problem that is solved using this handler.
Definition: BaseHandler.h:733
void copyContentsFromOtherHandler(const BaseHandler< BaseParticle > &BH)
Function that copies the contents (vector of pointers, maxObject_, nextId_, DPMBase_) from one handle...
Definition: BaseHandler.h:348
void clear() override
Empties the whole ParticleHandler by removing all BaseParticle.
Definition: ParticleHandler.cc:971
void computeSmallestParticle()
Computes the smallest particle (by interaction radius) and sets it in smallestParticle_.
Definition: ParticleHandler.cc:445
void computeLargestParticle()
Computes the largest particle (by interaction radius) and sets it in largestParticle_.
Definition: ParticleHandler.cc:469

References clear(), computeLargestParticle(), computeSmallestParticle(), BaseHandler< BaseParticle >::copyContentsFromOtherHandler(), DEBUG, BaseHandler< T >::getDPMBase(), largestParticle_, logger, BaseHandler< BaseParticle >::objects_, BaseHandler< BaseParticle >::setDPMBase(), and smallestParticle_.

◆ ~ParticleHandler()

ParticleHandler::~ParticleHandler ( )
override

Destructor, it destructs the ParticleHandler and all BaseParticle it contains.

Set the pointers to largestParticle_ and smallestParticle_ to nullptr, all BaseParticle are destroyed by the BaseHandler afterwards.

87 {
88  //First reset the pointers, such that they are not checked twice when removing particles
89  largestParticle_ = nullptr;
90  smallestParticle_ = nullptr;
91  logger(DEBUG, "ParticleHandler::~ParticleHandler() finished");
92 }

References DEBUG, largestParticle_, logger, and smallestParticle_.

Member Function Documentation

◆ actionsAfterTimeStep()

void ParticleHandler::actionsAfterTimeStep ( )
1379 {
1380  for (auto i: *this)
1381  {
1382  i->actionsAfterTimeStep();
1383  }
1384 }
int i
Definition: BiCGSTAB_step_by_step.cpp:9

References i.

Referenced by DPMBase::computeOneTimeStep().

◆ actionsBeforeTimeStep()

void ParticleHandler::actionsBeforeTimeStep ( )
1371 {
1372  for (auto i: *this)
1373  {
1374  i->actionsBeforeTimeStep();
1375  }
1376 }

References i.

Referenced by DPMBase::computeOneTimeStep().

◆ addedFixedParticle()

void ParticleHandler::addedFixedParticle ( )

Increment of the number of fixed particles.

Todo:
MX: For Jonny, is this still required, keeping the parallel code in mind?
1233 {
1234  NFixedParticles_++;
1235 }
unsigned int NFixedParticles_
Number of fixed particles.
Definition: ParticleHandler.h:466

References NFixedParticles_.

Referenced by BaseParticle::fixParticle().

◆ addExistingObject()

void ParticleHandler::addExistingObject ( BaseParticle P)
override

Adds a BaseParticle to the ParticleHandler.

Parameters
[in]PA pointer to the BaseParticle that has to be added.

To add a BaseParticle to the ParticleHandler, first check if it has a species, since it is as common bug to use a BaseParticle without species, which leads to a segmentation fault. To help the user with debugging, a warning is given if a particle without species is added. After that, the actions for adding the particle to the BaseHandler are taken, which include adding it to the vector of pointers to all BaseParticle and assigning the correct id and index. Then the particle is added to the HGrid, the particle is told that this is its handler, its mass is computed and finally it is checked if this is the smallest or largest particle in this ParticleHandler. The particle exists, in other words: it has been made before. This implies it already got an id Attached to it and hence we don't want to assign a new ID to it.

111 {
112  if (P->getSpecies() == nullptr)
113  {
114  logger(WARN, "WARNING: The particle with ID % that is added in ParticleHandler::addObject does not have a "
115  "species yet. Please make sure that you have "
116  "set the species somewhere in the driver code.", P->getId());
117  }
118 
119  //Puts the particle in the Particle list
121  if (getDPMBase() != nullptr)
122  {
123  //This places the particle in this grid
125  //This computes where the particle currently is in the grid
127  }
128  //set the particleHandler pointer
129  P->setHandler(this);
130  //compute mass of the particle
131  P->getSpecies()->computeMass(P);
132  //Check if this particle has new extrema
133  checkExtrema(P);
134 }
@ WARN
virtual void addExistingObject(T *O)
Adds an existing object to the BaseHandler without changing the id of the object.
Definition: BaseHandler.h:398
virtual void hGridUpdateParticle(BaseParticle *obj UNUSED)
Definition: DPMBase.cc:1693
virtual void hGridInsertParticle(BaseParticle *obj UNUSED)
Definition: DPMBase.cc:1686
void checkExtrema(BaseParticle *P)
Checks if the extrema of this ParticleHandler needs updating.
Definition: ParticleHandler.cc:1171
double P
Uniform pressure.
Definition: TwenteMeshGluing.cpp:77

References BaseHandler< T >::addExistingObject(), checkExtrema(), BaseHandler< BaseParticle >::getDPMBase(), DPMBase::hGridInsertParticle(), DPMBase::hGridUpdateParticle(), logger, Global_Physical_Variables::P, and WARN.

Referenced by readAndAddObject().

◆ addGhostObject() [1/2]

void ParticleHandler::addGhostObject ( BaseParticle P)
override

Adds a BaseParticle to the ParticleHandler.

This is a special function that is used in the parallel code. The functions differs from the standard parallel implementation as this function does not check if the particle is an mpi particle and the mpi domain is not updated. When the domain adds mpi particles to the simulation these checks are not required. It also doesnt increase the id of the particle

Parameters
[in]Apointer to the BaseParticle that has to be added.
337 {
338 #ifdef MERCURYDPM_USE_MPI
339  if (P->getSpecies() == nullptr)
340  {
341  logger(WARN, "[ParticleHandler::adGhostObject(BaseParticle*)] "
342  "WARNING: The particle with ID % that is added in "
343  "ParticleHandler::addObject does not have a species yet."
344  " Please make sure that you have "
345  "set the species somewhere in the driver code.", P->getId());
346  }
347 
348  MPIContainer& communicator = MPIContainer::Instance();
349  //Puts the particle in the Particle list
351  if (getDPMBase() != nullptr)
352  {
353  //This places the particle in this grid
355  //This computes where the particle currently is in the grid
357 
358  // broadcast particle addition
359  getDPMBase()->handleParticleAddition(P->getId(), P);
360  }
361  //set the particleHandler pointer
362  P->setHandler(this);
363  //compute mass of the particle
364  P->getSpecies()->computeMass(P) ;
365  //Check if this particle has new extrema
366  checkExtrema(P);
367 
368  if (getDPMBase() != nullptr)
369  {
370  // Broadcasts the existence of a new particle
371  getDPMBase()->handleParticleAddition(P->getId(), P);
372  }
373 #else
374  logger(INFO,
375  "Function ParticleHandler::mpiAddObject(BaseParticle* P) should only be called when compiling with parallel build on");
376 #endif
377 }
@ INFO
virtual void addGhostObject(T *O)
Adds a new Object to the BaseHandler. called by the to avoid increasing the id.
Definition: BaseHandler.h:425
virtual void handleParticleAddition(unsigned int id, BaseParticle *p)
Handles the addition of particles to the particleHandler.
Definition: DPMBase.cc:5566
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

References BaseHandler< T >::addGhostObject(), checkExtrema(), BaseHandler< BaseParticle >::getDPMBase(), DPMBase::handleParticleAddition(), DPMBase::hGridInsertParticle(), DPMBase::hGridUpdateParticle(), INFO, MPIContainer::Instance(), logger, Global_Physical_Variables::P, and WARN.

◆ addGhostObject() [2/2]

void ParticleHandler::addGhostObject ( int  fromProcessor,
int  toProcessor,
BaseParticle p 
)

Adds a ghost particle located at fromProcessor to toProcessor.

This function adds a ghost particle from one processor to another processor.

When a periodic ghost particle or a mpi ghost particle is created, the ID should not be updated and hence the addGhostObject needs to be called. Adding a ghost keeps the ID constant. Additionally this function copies the particle information given on processor fromProcessor to all other processors. The target processor then adds the particle.

Parameters
[in]fromProcessorprocessor containing the particle data
[in]toProcessorprocessor that needs to add the particle
[in,out]particlethat contains the data and receives the data
276 {
277 #ifdef MERCURYDPM_USE_MPI
278  MPIContainer& communicator = MPIContainer::Instance();
279  if (fromProcessor == toProcessor)
280  {
281  if (communicator.getProcessorID() == fromProcessor)
282  {
283  addGhostObject(p);
284  }
285 
286  return;
287  }
288 
289  //Create communication information object
291  MPIParticle pInfo;
292  int tag;
293  if (communicator.getProcessorID() == fromProcessor)
294  {
296  tag = fromProcessor*MAX_PROC*10 + toProcessor*10 + MercuryMPITag::PARTICLE_DATA;
297  communicator.send(&pInfo, MercuryMPIType::PARTICLE, 1, toProcessor, tag);
298  }
299 
300  if (communicator.getProcessorID() == toProcessor)
301  {
302  tag = fromProcessor*MAX_PROC*10 + toProcessor*10 + MercuryMPITag::PARTICLE_DATA;
303  communicator.receive(&pInfo, MercuryMPIType::PARTICLE, 1, fromProcessor, tag);
304  }
305 
306  //Sync the communication
307  communicator.sync();
308 
309  //Add the data to the new particle
310  if (communicator.getProcessorID() == toProcessor)
311  {
312  copyDataFromMPIParticleToParticle(&pInfo, p, this);
313  }
314 
315 
316  //Only toProcessor adds the particle, quietly without any check with other processors
317  //IMPORTANT NOTE: When adding a periodic particle in parallel, this is performed just before
318  //finding new particles. For that reason we dont have to add a ghostparticle to the mpi communication lists
319  if (communicator.getProcessorID() == toProcessor)
320  {
321  addGhostObject(p);
322  }
323 #else
324  logger(WARN,
325  "Function addGhostObject(int fromProcessor, int toProcessor, BaseParticle* p) should not be used in serial code");
326 #endif
327 }
#define MAX_PROC
Definition: GeneralDefine.h:30
@ PARTICLE_DATA
Definition: MpiContainer.h:58
@ PARTICLE
Definition: MpiContainer.h:46
void copyDataFromMPIParticleToParticle(MPIParticle *bP, BaseParticle *p, ParticleHandler *particleHandler)
Copies data from an MPIParticle class to a BaseParticle and sets the particleHandler and species.
Definition: MpiDataClass.cc:84
float * p
Definition: Tutorial_Map_using.cpp:9
std::enable_if< std::is_scalar< T >::value, void >::type receive(T &t, int from, int tag)
asynchronously receive a scalar from some other processor.
Definition: MpiContainer.h:200
void sync()
Process all pending asynchronous communication requests before continuing.
Definition: MpiContainer.h:131
std::size_t getProcessorID()
Reduces a scalar on all processors to one scalar on a target processor.
Definition: MpiContainer.cc:92
std::enable_if< std::is_scalar< T >::value, void >::type send(T &t, int to, int tag)
Asynchronously send a scalar to some other processor.
Definition: MpiContainer.h:150
Data class to send a particle over MPI.
Definition: MpiDataClass.h:60
void copyDataFromParticleToMPIParticle(BaseParticle *p)
Definition: MpiDataClass.cc:110
void addGhostObject(int fromProcessor, int toProcessor, BaseParticle *p)
Adds a ghost particle located at fromProcessor to toProcessor.
Definition: ParticleHandler.cc:275

References copyDataFromMPIParticleToParticle(), MPISphericalParticle::copyDataFromParticleToMPIParticle(), MPIContainer::getProcessorID(), MPIContainer::Instance(), logger, MAX_PROC, p, PARTICLE, PARTICLE_DATA, MPIContainer::receive(), MPIContainer::send(), MPIContainer::sync(), and WARN.

Referenced by PeriodicBoundaryHandler::processLocalGhostParticles(), Domain::processReceivedBoundaryParticleData(), and PeriodicBoundaryHandler::processReceivedGhostParticleData().

◆ addObject() [1/2]

void ParticleHandler::addObject ( BaseParticle P)
override

Adds a BaseParticle to the ParticleHandler.

Parameters
[in]PA pointer to the BaseParticle that has to be added.

To add a BaseParticle to the ParticleHandler, first check if it has a species, since it is as common bug to use a BaseParticle without species, which leads to a segmentation fault. To help the user with debugging, a warning is given if a particle without species is added. After that, the actions for adding the particle to the BaseHandler are taken, which include adding it to the vector of pointers to all BaseParticle and assigning the correct id and index. Then the particle is added to the HGrid, the particle is told that this is its handler, its mass is computed and finally it is checked if this is the smallest or largest particle in this ParticleHandler.

151 {
152  if (P->getSpecies() == nullptr)
153  {
154  logger(WARN, "WARNING: The particle with ID % that is added in "
155  "ParticleHandler::addObject does not have a species yet. "
156  "Please make sure that you have "
157  "set the species somewhere in the driver code.", P->getId());
158  }
159 #ifdef MERCURYDPM_USE_MPI
160  bool insertParticle;
161  //Check if the particle P should be added to the current domain
162  if (NUMBER_OF_PROCESSORS == 1)
163  {
164  insertParticle = true;
165  }
166  else
167  {
168  insertParticle = getDPMBase()->mpiInsertParticleCheck(P);
169  }
170 
171  //Add the particle if it is in the mpi domain or if the domain is not yet defined
172  if(insertParticle)
173  {
174 #endif
175  //Puts the particle in the Particle list
177  if (getDPMBase() != nullptr)
178  {
179  //This places the particle in this grid
181  //This computes where the particle currently is in the grid
183 
184  // Broadcasts the existance of a new particle
185  getDPMBase()->handleParticleAddition(P->getId(), P);
186 
187  // Store the time stamp of creation.
188  P->setTimeStamp(getDPMBase()->getNumberOfTimeSteps());
189  }
190  //set the particleHandler pointer
191  P->setHandler(this);
192  //compute mass of the particle
193  P->getSpecies()->computeMass(P);
194  //Check if this particle has new extrema
195  checkExtrema(P);
196  if (!P->isSphericalParticle())
197  {
198  getDPMBase()->setRotation(true);
199  }
200 
201  P->actionsAfterAddObject();
202 
203  if (getDPMBase() != nullptr)
204  {
205  // Broadcasts the existence of a new particle
206  getDPMBase()->handleParticleAddition(P->getId(), P);
207  }
208 
209 #ifdef MERCURYDPM_USE_MPI
210  P->setPeriodicComplexity(std::vector<int>(0));
211  }
212  else
213  {
214  logger.assert_debug(!P->isMPIParticle(),"Can't add mpi particle as it does not exist");
215  logger.assert_debug(!P->isPeriodicGhostParticle(),"Can't add mpi particle as it does not exist");
216  //Somehwere a really new particle has been added, so to keep the ID's globally unique, we also update
217  //the unique value on all other processors
219  }
220 
221  //Add the particle to the ghost particle lists
223  //Check if this new particle requires an update in the mpi grid (interactionDistance).
225 
226  //Delete the particle that was supposed to be added (but was not)
227  if(!insertParticle)
228  {
229  delete P;
230  }
231 #endif
232 }
#define NUMBER_OF_PROCESSORS
For the MPI communication routines this quantity is often required. defining this macro makes the cod...
Definition: GeneralDefine.h:41
virtual void addObject(T *object)
Adds a new Object to the BaseHandler.
Definition: BaseHandler.h:412
void increaseId()
Definition: BaseHandler.h:232
void setRotation(bool rotation)
Sets whether particle rotation is enabled or disabled.
Definition: DPMBase.h:547
void insertGhostParticle(BaseParticle *P)
This function inserts a particle in the mpi communication boundaries.
Definition: DPMBase.cc:1800
void updateGhostGrid(BaseParticle *P)
Checks if the Domain/periodic interaction distance needs to be updated and updates it accordingly.
Definition: DPMBase.cc:1826
bool mpiInsertParticleCheck(BaseParticle *P)
Function that checks if the mpi particle should really be inserted by the current domain.
Definition: DPMBase.cc:1721

References BaseHandler< T >::addObject(), checkExtrema(), BaseHandler< BaseParticle >::getDPMBase(), DPMBase::handleParticleAddition(), DPMBase::hGridInsertParticle(), DPMBase::hGridUpdateParticle(), BaseHandler< T >::increaseId(), DPMBase::insertGhostParticle(), logger, DPMBase::mpiInsertParticleCheck(), NUMBER_OF_PROCESSORS, Global_Physical_Variables::P, DPMBase::setRotation(), DPMBase::updateGhostGrid(), and WARN.

Referenced by ChuteWithPeriodicInflow::AddContinuingBottom(), addObject(), CircularPeriodicBoundary::checkBoundaryAfterParticleMoved(), ConstantMassFlowMaserBoundary::checkBoundaryAfterParticleMoved(), SubcriticalMaserBoundary::checkBoundaryAfterParticleMoved(), SubcriticalMaserBoundaryTEST::checkBoundaryAfterParticleMoved(), ChuteWithContraction::ChuteWithContraction(), ChuteWithPeriodicInflowAndContinuingBottom::ChuteWithPeriodicInflowAndContinuingBottom(), ChuteWithPeriodicInflowAndContraction::ChuteWithPeriodicInflowAndContraction(), ChuteWithPeriodicInflowAndVariableBottom::ChuteWithPeriodicInflowAndVariableBottom(), ContractionWithPeriodicInflow::ContractionWithPeriodicInflow(), SubcriticalMaserBoundaryTEST::copyExtraParticles(), CurvyChute::createBottom(), PeriodicBoundary::createGhostParticle(), TimeDependentPeriodicBoundary::createGhostParticle(), LeesEdwardsBoundary::createHorizontalPeriodicParticle(), ShearBoxBoundary::createHorizontalPeriodicParticle(), AngledPeriodicBoundary::createPeriodicParticle(), CircularPeriodicBoundary::createPeriodicParticle(), ConstantMassFlowMaserBoundary::createPeriodicParticle(), SubcriticalMaserBoundary::createPeriodicParticle(), LeesEdwardsBoundary::createVerticalPeriodicParticle(), ShearBoxBoundary::createVerticalPeriodicParticle(), SubcriticalMaserBoundaryTEST::extendBottom(), ChuteWithPeriodicInflow::ExtendInWidth(), ChuteWithPeriodicInflow::integrateBeforeForceComputation(), and ContactDetectionTester::setupParticles().

◆ addObject() [2/2]

void ParticleHandler::addObject ( int  fromProcessor,
BaseParticle p 
)

Adds a BaseParticle located at processor fromProcessor to toProcessor.

This function adds a particle to the simulation where the information of the particle is not available by the target processor.

When a certain processor generates a particle which needs to be inserted in a domain of another processor, this information needs to be transfered before the particle can be actually added

Parameters
[in]fromProcessorprocessor containing the particle data
[in,out]particlethat contains the data and receives the data
243 {
244 #ifdef MERCURYDPM_USE_MPI
245  MPIContainer& communicator = MPIContainer::Instance();
246 
247  //The processor that contains the particle that needs to be copied needs to identify the target, and communicate this
248  MPIParticle pInfo;
249  if (communicator.getProcessorID() == fromProcessor)
250  {
252  }
253 
254  //Broadcast from processor i
255  communicator.broadcast(&pInfo,MercuryMPIType::PARTICLE,fromProcessor);
256  copyDataFromMPIParticleToParticle(&pInfo, p, this);
257 
258  //All processors now know have the same information and we can proceed with the collective add
259  addObject(p);
260 #else
261  logger(WARN, "Function addObject(int fromProcessor, BaseParticle* p) should not be used in serial code");
262 #endif
263 }
std::enable_if< std::is_scalar< T >::value, void >::type broadcast(T &t, int fromProcessor=0)
Broadcasts a scalar from the root to all other processors.
Definition: MpiContainer.h:420
void addObject(BaseParticle *P) override
Adds a BaseParticle to the ParticleHandler.
Definition: ParticleHandler.cc:150

References addObject(), MPIContainer::broadcast(), copyDataFromMPIParticleToParticle(), MPISphericalParticle::copyDataFromParticleToMPIParticle(), MPIContainer::getProcessorID(), MPIContainer::Instance(), logger, p, PARTICLE, and WARN.

◆ checkExtrema()

void ParticleHandler::checkExtrema ( BaseParticle P)

Checks if the extrema of this ParticleHandler needs updating.

 \param[in] type The first value of the position.
 \param[in] is The input stream from which the information is read.
 \details The old objects did not have their type in the beginning of the line.
          Instead, the first string of the file was the position in x-direction.
          Since we already read the first string of the file, we need to give
          it to this function and convert it to the position in x-direction.
          The rest of the stream is then read in the usual way.
&zwj;/

void ParticleHandler::readAndCreateOldObject(std::istream& is, const std::string& type) { //read in next line std::stringstream line; helpers::getLineFromStringStream(is, line); logger(VERBOSE, line.str()); //stdcout << line.str() << std::endl;

BaseParticle particle;

//Declare all properties of the particle unsigned int indSpecies; Mdouble radius, inverseMass, inverseInertia; Vec3D position, velocity, euler, angularVelocity;

//Read all values position.X = atof(type.c_str());

line >> position.Y >> position.Z >> velocity >> radius >> euler >> angularVelocity >> inverseMass >> inverseInertia >> indSpecies;

//Put the values in the particle particle.setSpecies(getDPMBase()->speciesHandler.getObject(indSpecies)); particle.setPosition(position); particle.setVelocity(velocity); particle.setRadius(radius); Quaternion q; q.setEuler(euler); particle.setOrientation(q); particle.setAngularVelocity(angularVelocity); if (inverseMass == 0.0) particle.fixParticle(); else { particle.setInverseInertia(MatrixSymmetric3D(1,0,0,1,0,1)*inverseInertia); }

//Put the particle in the handler copyAndAddObject(particle); }

void ParticleHandler::write(std::ostream& os) const {

os << "Particles " << getSize() << '\n';
for (BaseParticle* it : *this)
{
    os << (*it) << '\n';
}

os.flush(); }

/*!

Parameters
[in]PA pointer to the particle, which properties have to be checked against the ParticleHandlers extrema.
1172 {
1173  if (P == largestParticle_)
1174  {
1175  //if the properties of the largest particle changes
1177  }
1178  else if (!largestParticle_ || P->getMaxInteractionRadius() > largestParticle_->getMaxInteractionRadius())
1179  {
1180  largestParticle_ = P;
1181  }
1182 
1183  if (P == smallestParticle_)
1184  {
1185  //if the properties of the smallest particle changes
1187  }
1188  else if (!smallestParticle_ || P->getMaxInteractionRadius() < smallestParticle_->getMaxInteractionRadius())
1189  {
1190  smallestParticle_ = P;
1191  }
1192 }
Mdouble getMaxInteractionRadius() const
Returns the particle's interaction radius, which might be different from radius_ (e....
Definition: BaseParticle.h:345

References computeLargestParticle(), computeSmallestParticle(), BaseParticle::getMaxInteractionRadius(), largestParticle_, Global_Physical_Variables::P, and smallestParticle_.

Referenced by addExistingObject(), addGhostObject(), addObject(), and BaseParticle::setRadius().

◆ checkExtremaOnDelete()

void ParticleHandler::checkExtremaOnDelete ( BaseParticle P)

Checks if the extrema of this ParticleHandler needs updating when a particle is deleted.

Parameters
[in]PA pointer to the particle, which is going to get deleted.
1198 {
1199  if (P == largestParticle_)
1200  {
1202  }
1203  if (P == smallestParticle_)
1204  {
1206  }
1207 }

References computeLargestParticle(), computeSmallestParticle(), largestParticle_, Global_Physical_Variables::P, and smallestParticle_.

Referenced by BaseParticle::~BaseParticle(), and SuperQuadricParticle::~SuperQuadricParticle().

◆ clear()

void ParticleHandler::clear ( )
overridevirtual

Empties the whole ParticleHandler by removing all BaseParticle.

Note that the pointers to smallestParticle_ and largestParticle_ are set to nullptr since these particles don't exist anymore after calling this function.

Reimplemented from BaseHandler< BaseParticle >.

972 {
973  smallestParticle_ = nullptr;
974  largestParticle_ = nullptr;
975 
976  for (auto p0: *this)
977  {
978  getDPMBase()->handleParticleRemoval(p0->getId());
979  }
980 
982 }
Vector3f p0
Definition: MatrixBase_all.cpp:2
virtual void clear()
Empties the whole BaseHandler by removing all Objects and setting all other variables to 0.
Definition: BaseHandler.h:536
virtual void handleParticleRemoval(unsigned int id)
Handles the removal of particles from the particleHandler.
Definition: DPMBase.cc:5551

References BaseHandler< T >::clear(), BaseHandler< BaseParticle >::getDPMBase(), DPMBase::handleParticleRemoval(), largestParticle_, p0, and smallestParticle_.

Referenced by NautaMixer::addParticles(), NautaMixer::addParticlesAtWall(), CGExactOverlapUnitTest::clear(), CGHandler::evaluateRestartFiles(), HorizontalMixer::introduceParticlesAtWall(), HorizontalMixer::introduceParticlesInDomain(), HorizontalMixer::introduceSingleParticle(), operator=(), ParticleHandler(), FileReader::read(), DPMBase::read(), DPMBase::readNextDataFile(), HorizontalMixerWalls::setupInitialConditions(), FreeCooling2DinWallsDemo::setupInitialConditions(), FreeCooling3DDemoProblem::setupInitialConditions(), FreeCooling3DinWallsDemo::setupInitialConditions(), FreeCoolingDemoProblem::setupInitialConditions(), CoilSelfTest::setupInitialConditions(), ParticleParticleInteraction::setupInitialConditions(), ParticleParticleInteractionWithPlasticForces::setupInitialConditions(), ParticleWallInteraction::setupInitialConditions(), EnergyUnitTest::setupInitialConditions(), HertzianSinterForceUnitTest::setupInitialConditions(), MD_demo::setupInitialConditions(), PlasticForceUnitTest::setupInitialConditions(), SeparateFilesSelfTest::setupInitialConditions(), SinterForceUnitTest::setupInitialConditions(), and BaseCluster::setupInitialConditions().

◆ computeAllMasses() [1/2]

void ParticleHandler::computeAllMasses ( )

Computes the mass for all BaseParticle in this ParticleHandler.

1225 {
1226  for (BaseParticle* particle : objects_)
1227  {
1228  particle->getSpecies()->computeMass(particle);
1229  }
1230 }
Definition: BaseParticle.h:33

References BaseHandler< BaseParticle >::objects_.

◆ computeAllMasses() [2/2]

void ParticleHandler::computeAllMasses ( unsigned int  indSpecies)

Computes the mass for all BaseParticle of the given species in this ParticleHandler.

Parameters
[in]indSpeciesUnsigned integer with the index of the species for which the masses must be computed.
1214 {
1215  for (BaseParticle* particle : objects_)
1216  {
1217  if (particle->getIndSpecies() == indSpecies)
1218  {
1219  particle->getSpecies()->computeMass(particle);
1220  }
1221  }
1222 }

References BaseHandler< BaseParticle >::objects_.

Referenced by SpeciesHandler::addObject(), DPMBase::initialiseSolve(), DPMBase::readNextDataFile(), ParticleSpecies::setDensity(), and DPMBase::setParticleDimensions().

◆ computeLargestParticle()

void ParticleHandler::computeLargestParticle ( )

Computes the largest particle (by interaction radius) and sets it in largestParticle_.

470 {
471  if (getSize() == 0)
472  {
473  logger(DEBUG, "No particles, so cannot compute the largest particle.");
474  largestParticle_ = nullptr;
475  return;
476  }
478  largestParticle_ = nullptr;
479  for (BaseParticle* const particle : objects_)
480  {
481  //if (!(particle->isMPIParticle() || particle->isPeriodicGhostParticle()))
482  {
483  if (particle->getMaxInteractionRadius() > max)
484  {
485  max = particle->getMaxInteractionRadius();
486  largestParticle_ = particle;
487  }
488  }
489  }
490 }
unsigned int getSize() const
Gets the size of the particleHandler (including mpi and periodic particles)
Definition: BaseHandler.h:663
#define max(a, b)
Definition: datatypes.h:23

References DEBUG, BaseHandler< BaseParticle >::getSize(), largestParticle_, logger, max, and BaseHandler< BaseParticle >::objects_.

Referenced by DeletionBoundary::checkBoundaryAfterParticleMoved(), checkExtrema(), checkExtremaOnDelete(), operator=(), and ParticleHandler().

◆ computeSmallestParticle()

void ParticleHandler::computeSmallestParticle ( )

Computes the smallest particle (by interaction radius) and sets it in smallestParticle_.

446 {
447  if (getSize() == 0)
448  {
449  logger(DEBUG, "No particles, so cannot compute the smallest particle.");
450 
451  smallestParticle_ = nullptr;
452  return;
453  }
455  smallestParticle_ = nullptr;
456  for (BaseParticle* const particle : objects_)
457  {
458  //if (!(particle->isMPIParticle() || particle->isPeriodicGhostParticle()))
459  {
460  if (particle->getMaxInteractionRadius() < min)
461  {
462  min = particle->getMaxInteractionRadius();
463  smallestParticle_ = particle;
464  }
465  }
466  }
467 }
#define min(a, b)
Definition: datatypes.h:22

References DEBUG, BaseHandler< BaseParticle >::getSize(), logger, max, min, BaseHandler< BaseParticle >::objects_, and smallestParticle_.

Referenced by DeletionBoundary::checkBoundaryAfterParticleMoved(), checkExtrema(), checkExtremaOnDelete(), operator=(), and ParticleHandler().

◆ createObject()

BaseParticle * ParticleHandler::createObject ( const std::string &  type)
static

Reads BaseParticle into the ParticleHandler from restart data.

Parameters
[in]isThe input stream from which the information is read.
1020 {
1021  if (type == "BaseParticle") {
1022  //for backwards compatibility
1023  return new SphericalParticle;
1024  }
1025  else if (type == "SphericalParticle")
1026  {
1027  return new SphericalParticle;
1028  }
1029  else if (type == "LiquidFilmBaseParticle" || type == "LiquidFilmParticle")
1030  {
1031  return new LiquidFilmParticle;
1032  }
1033  else if (type == "SuperQuadricParticle")
1034  {
1035  return new SuperQuadricParticle;
1036  }
1037  else if (type == "ThermalParticle")
1038  {
1039  return new ThermalParticle;
1040  }
1041  else if (type == "MeltableParticle")
1042  {
1043  return new MeltableParticle;
1044  }
1045  else if (type == "HeatFluidCoupledBaseParticle")
1046  {
1047  return new HeatFluidCoupledParticle;
1048  }
1049  else if (type == "Clump")
1050  {
1051  return new ClumpParticle;
1052  }
1053  else
1054  {
1055  logger(WARN, "Particle type % not understood in restart file. Particle will not be read.", type);
1056  return nullptr;
1057  }
1058 }
HeatFluidCoupled< SphericalParticle > HeatFluidCoupledParticle
Template specialisation of HeatFluidCoupled<Particle> for spherical particles.
Definition: HeatFluidCoupledParticle.h:111
LiquidFilm< SphericalParticle > LiquidFilmParticle
Definition: LiquidMigrationLSInteraction.h:22
Thermal< SphericalParticle > ThermalParticle
Definition: ThermalParticle.h:165
Definition: ClumpParticle.h:20
Definition: MeltableParticle.h:15
A spherical particle is the most simple particle used in MercuryDPM.
Definition: SphericalParticle.h:16
Definition: SuperQuadricParticle.h:36
type
Definition: compute_granudrum_aor.py:141

References logger, compute_granudrum_aor::type, and WARN.

Referenced by readAndCreateObject().

◆ getAngularMomentum()

Vec3D ParticleHandler::getAngularMomentum ( ) const
670 {
671  Vec3D momentum = {0, 0, 0};
672  for (auto p : *this)
673  if (!(p->isFixed() || p->isMPIParticle() || p->isPeriodicGhostParticle()))
674  momentum += p->getAngularMomentum();
675  return getMPISum(momentum);
676 }
Vec3D getMPISum(Vec3D &val)
Definition: MpiDataClass.cc:178
Definition: Kernel/Math/Vector.h:30

References getMPISum(), and p.

◆ getCentreOfMass()

Vec3D ParticleHandler::getCentreOfMass ( ) const
649 {
650  Mdouble m = getMass();
651  if (m == 0)
652  {
654  return nanvec;
655  }
656  else
657  return getMassTimesPosition() / m;
658 }
Vec3D getMassTimesPosition() const
Definition: ParticleHandler.cc:630
Mdouble getMass() const
Definition: ParticleHandler.cc:605
int * m
Definition: level2_cplx_impl.h:294
const Mdouble NaN
Definition: GeneralDefine.h:22

References getMass(), getMassTimesPosition(), m, and constants::NaN.

Referenced by DPMBase::getCentreOfMass(), and RotatingDrumWet::printTime().

◆ getFastestParticle()

BaseParticle * ParticleHandler::getFastestParticle ( ) const
705 {
706 #ifdef MERCURYDPM_USE_MPI
707  logger(ERROR,"This function should not be used in parallel");
708 #endif
709  return getFastestParticleLocal();
710 }
@ ERROR
BaseParticle * getFastestParticleLocal() const
Gets a pointer to the fastest BaseParticle in this ParticleHandler.
Definition: ParticleHandler.cc:681

References ERROR, getFastestParticleLocal(), and logger.

◆ getFastestParticleLocal()

BaseParticle * ParticleHandler::getFastestParticleLocal ( ) const

Gets a pointer to the fastest BaseParticle in this ParticleHandler.

Returns
A pointer to the fastest BaseParticle in this ParticleHandler.
682 {
683  if (getSize() == 0)
684  {
685  logger(WARN, "No particles to set getFastestParticle()");
686  return nullptr;
687  }
688  BaseParticle* p = nullptr;
690  for (BaseParticle* const pLoop : objects_)
691  {
692  if (!(pLoop->isMPIParticle() || pLoop->isPeriodicGhostParticle()))
693  {
694  if ((pLoop->getVelocity().getLength()) > maxSpeed)
695  {
696  maxSpeed = pLoop->getVelocity().getLength();
697  p = pLoop;
698  }
699  }
700  }
701  return p;
702 }

References BaseHandler< BaseParticle >::getSize(), logger, max, BaseHandler< BaseParticle >::objects_, p, and WARN.

Referenced by getFastestParticle().

◆ getHighestPositionComponentParticle()

BaseParticle * ParticleHandler::getHighestPositionComponentParticle ( int  i) const

Gets a pointer to the particle with the highest coordinates in direction i in this ParticleHandler.

887 {
888 #ifdef MERCURYDPM_USE_MPI
889  logger(ERROR,"This function should not be used in parallel");
890 #endif
892 }
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

References ERROR, getHighestPositionComponentParticleLocal(), i, and logger.

◆ getHighestPositionComponentParticleLocal()

BaseParticle * ParticleHandler::getHighestPositionComponentParticleLocal ( int  i) const

Gets a pointer to the particle with the highest coordinates in direction i in this ParticleHandler.

Parameters
[in]iDirection for which one wants the particle with highest coordinates.
Returns
A pointer to the particle with the highest coordinates in direction i in this ParticleHandler.
863 {
864  if (getSize() == 0)
865  {
866  logger(WARN, "No getHighestPositionComponentParticle(const int i) since there are no particles.");
867  return nullptr;
868  }
869  BaseParticle* p = nullptr;
871  for (BaseParticle* const pLoop : objects_)
872  {
873  if (!(pLoop->isMPIParticle() || pLoop->isPeriodicGhostParticle()))
874  {
875  if (pLoop->getPosition().getComponent(i) > max)
876  {
877  max = pLoop->getPosition().getComponent(i);
878  p = pLoop;
879  }
880  }
881  }
882 
883  return p;
884 }

References BaseHandler< BaseParticle >::getSize(), i, logger, max, BaseHandler< BaseParticle >::objects_, p, and WARN.

Referenced by getHighestPositionComponentParticle().

◆ getHighestPositionX()

Mdouble ParticleHandler::getHighestPositionX ( ) const

Function returns the highest position in the x-direction.

This is a prototype example of how to obtain global particle attributes in the parallel code

Returns
Returns the highest x-position value of all particles
990 {
991  //Define the attribute function
992  std::function<Mdouble(BaseParticle*)> particleAttribute = [](BaseParticle* p) { return p->getPosition().X; };
993 
994  //Obtain the MAX attribute
995  Mdouble positionXGlobal = getParticleAttribute<Mdouble>(particleAttribute, AttributeType::MAX);
996 
997  return positionXGlobal;
998 }
double Mdouble
Definition: GeneralDefine.h:13
@ MAX
Definition: ParticleHandler.h:18

References MAX, and p.

◆ getHighestVelocityComponentParticle()

BaseParticle * ParticleHandler::getHighestVelocityComponentParticle ( int  i) const

Gets a pointer to the particle with the highest velocity in direction i in this ParticleHandler.

959 {
960 #ifdef MERCURYDPM_USE_MPI
961  logger(ERROR,"This function should not be used in parallel");
962 #endif
964 }
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

References ERROR, getHighestVelocityComponentParticleLocal(), i, and logger.

◆ getHighestVelocityComponentParticleLocal()

BaseParticle * ParticleHandler::getHighestVelocityComponentParticleLocal ( int  i) const

Gets a pointer to the particle with the highest velocity in direction i in this ParticleHandler.

Parameters
[in]iDirection for which you want the particle with highest velocity.
Returns
A pointer to the particle with the highest velocity in direction i in this ParticleHandler.
936 {
937  if (!getSize())
938  {
939  logger(WARN, "No getHighestVelocityComponentParticle(const int i) since there are no particles");
940  return nullptr;
941  }
942  BaseParticle* p = nullptr;
944  for (BaseParticle* const pLoop : objects_)
945  {
946  if (!(pLoop->isMPIParticle() || pLoop->isPeriodicGhostParticle()))
947  {
948  if (pLoop->getVelocity().getComponent(i) > max)
949  {
950  max = pLoop->getVelocity().getComponent(i);
951  p = pLoop;
952  }
953  }
954  }
955  return p;
956 }

References BaseHandler< BaseParticle >::getSize(), i, logger, max, BaseHandler< BaseParticle >::objects_, p, and WARN.

Referenced by getHighestVelocityComponentParticle().

◆ getKineticEnergy()

Mdouble ParticleHandler::getKineticEnergy ( ) const
552 {
553 #ifdef MERCURYDPM_USE_MPI
554  Mdouble kineticEnergyLocal = getKineticEnergyLocal();
555  Mdouble kineticEnergyGlobal = 0.0;
556 
557  //sum up over all domains
558  MPIContainer& communicator = MPIContainer::Instance();
559  communicator.allReduce(kineticEnergyLocal, kineticEnergyGlobal, MPI_SUM);
560 
561  return kineticEnergyGlobal;
562 #else
563  return getKineticEnergyLocal();
564 #endif
565 }
Mdouble getKineticEnergyLocal() const
Definition: ParticleHandler.cc:538

References getKineticEnergyLocal(), and MPIContainer::Instance().

Referenced by DPMBase::writeEneTimeStep(), and LawinenBox::writeEneTimeStep().

◆ getKineticEnergyLocal()

Mdouble ParticleHandler::getKineticEnergyLocal ( ) const
private
539 {
540  Mdouble ene = 0;
541  for (auto p : *this)
542  {
543  if (!(p->isMPIParticle() || p->isPeriodicGhostParticle()))
544  {
545  ene += p->getKineticEnergy();
546  }
547  }
548  return ene;
549 }

References p.

Referenced by getKineticEnergy().

◆ getLargestInteractionRadius()

Mdouble ParticleHandler::getLargestInteractionRadius ( ) const

Returns the largest interaction radius.

767 {
768 #ifdef MERCURYDPM_USE_MPI
769  Mdouble largestInteractionRadiusLocal = getLargestInteractionRadiusLocal();
770  Mdouble largestInteractionRadiusGlobal = 0.0;
771 
772  //Obtain the global value
773  MPIContainer& communicator = MPIContainer::Instance();
774  communicator.allReduce(largestInteractionRadiusLocal, largestInteractionRadiusGlobal, MPI_MAX);
775 
776  return largestInteractionRadiusGlobal;
777 #else
779 #endif
780 }
Mdouble getLargestInteractionRadiusLocal() const
Returns the largest interaction radius of the current domain.
Definition: ParticleHandler.cc:751

References getLargestInteractionRadiusLocal(), and MPIContainer::Instance().

Referenced by Contact::actionsAfterSolve(), getPSD(), and DPMBase::read().

◆ getLargestInteractionRadiusLocal()

Mdouble ParticleHandler::getLargestInteractionRadiusLocal ( ) const

Returns the largest interaction radius of the current domain.

752 {
753  if (!(getLargestParticleLocal() == nullptr))
754  {
756  }
757  else
758  {
759  return 0.0;
760  }
761 }
BaseParticle * getLargestParticleLocal() const
Gets a pointer to the largest BaseParticle (by interactionRadius) in the ParticleHandler of the local...
Definition: ParticleHandler.cc:520
BaseParticle * getLargestParticle() const
Returns the pointer of the largest particle in the particle handler. When mercury is running in paral...
Definition: ParticleHandler.cc:528

References getLargestParticle(), getLargestParticleLocal(), and BaseParticle::getMaxInteractionRadius().

Referenced by MercuryBase::getHGridTargetMaxInteractionRadius(), and getLargestInteractionRadius().

◆ getLargestParticle()

BaseParticle * ParticleHandler::getLargestParticle ( ) const

◆ getLargestParticleLocal()

BaseParticle * ParticleHandler::getLargestParticleLocal ( ) const

Gets a pointer to the largest BaseParticle (by interactionRadius) in the ParticleHandler of the local domain.

Returns
A pointer to the largest BaseParticle (by interactionRadius) in this ParticleHandler.
521 {
522  return largestParticle_;
523 }

References largestParticle_.

Referenced by ConstantMassFlowMaserBoundary::activateMaser(), ConstantMassFlowMaserBoundary::createPeriodicParticle(), getLargestInteractionRadiusLocal(), getLargestParticle(), RotatingDrum::setupInitialConditions(), and DPMBase::writeFstatHeader().

◆ getLiquidFilmVolume()

double ParticleHandler::getLiquidFilmVolume ( ) const
1387 {
1388  double liquidVolume = 0;
1389  for (auto i : objects_) {
1390  auto j = dynamic_cast<LiquidFilmParticle*>(i);
1391  if (j and !j->isMPIParticle()) liquidVolume += j->getLiquidVolume();
1392  }
1393  return getMPISum(liquidVolume);
1394 };
Definition: LiquidFilmParticle.h:15
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References getMPISum(), i, j, and BaseHandler< BaseParticle >::objects_.

Referenced by LiquidMigrationMPI2Test::printTime().

◆ getLowestPositionComponentParticle()

BaseParticle * ParticleHandler::getLowestPositionComponentParticle ( int  i) const

Gets a pointer to the particle with the lowest coordinates in direction i in this ParticleHandler.

850 {
851 #ifdef MERCURYDPM_USE_MPI
852  logger(ERROR,"This function should not be used in parallel");
853 #endif
855 }
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

References ERROR, getLowestPositionComponentParticleLocal(), i, and logger.

◆ getLowestPositionComponentParticleLocal()

BaseParticle * ParticleHandler::getLowestPositionComponentParticleLocal ( int  i) const

Gets a pointer to the particle with the lowest coordinates in direction i in this ParticleHandler.

Parameters
[in]iDirection for which you want the particle with lowest coordinates.
Returns
A pointer to the particle with the lowest coordinates in the given direction in this ParticleHandler.
824 {
825 #ifdef MERCURYDPM_USE_MPI
826  logger(ERROR,"getLowestPositionComponentParticle() not implemented yet in parallel");
827 #endif
828  if (getSize() == 0)
829  {
830  logger(WARN, "No getLowestPositionComponentParticle(const int i) since there are no particles.");
831  return nullptr;
832  }
833  BaseParticle* p = nullptr;
835  for (BaseParticle* const pLoop : objects_)
836  {
837  if (!(pLoop->isMPIParticle() || pLoop->isPeriodicGhostParticle()))
838  {
839  if (pLoop->getPosition().getComponent(i) < min)
840  {
841  min = pLoop->getPosition().getComponent(i);
842  p = pLoop;
843  }
844  }
845  }
846  return p;
847 }

References ERROR, BaseHandler< BaseParticle >::getSize(), i, logger, max, min, BaseHandler< BaseParticle >::objects_, p, and WARN.

Referenced by getLowestPositionComponentParticle().

◆ getLowestVelocityComponentParticle()

BaseParticle * ParticleHandler::getLowestVelocityComponentParticle ( int  i) const

Gets a pointer to the particle with the lowest velocity in direction i in this ParticleHandler.

923 {
924 #ifdef MERCURYDPM_USE_MPI
925  logger(ERROR,"This function should not be used in parallel");
926 #endif
928 }
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

References ERROR, getLowestVelocityComponentParticleLocal(), i, and logger.

◆ getLowestVelocityComponentParticleLocal()

BaseParticle * ParticleHandler::getLowestVelocityComponentParticleLocal ( int  i) const

Gets a pointer to the particle with the lowest velocity in direction i in this ParticleHandler.

Parameters
[in]iDirection for which you want the particle with lowest velocity.
Returns
A pointer to the particle with the lowest velocity in direction i in this ParticleHandler.
900 {
901  if (getSize() == 0)
902  {
903  logger(WARN, "No getLowestVelocityComponentParticle(const int i) since there are no particles");
904  return nullptr;
905  }
906  BaseParticle* p = nullptr;
908  for (BaseParticle* const pLoop : objects_)
909  {
910  if (!(pLoop->isMPIParticle() || pLoop->isPeriodicGhostParticle()))
911  {
912  if (pLoop->getVelocity().getComponent(i) < min)
913  {
914  min = pLoop->getVelocity().getComponent(i);
915  p = pLoop;
916  }
917  }
918  }
919  return p;
920 }

References BaseHandler< BaseParticle >::getSize(), i, logger, max, min, BaseHandler< BaseParticle >::objects_, p, and WARN.

Referenced by getLowestVelocityComponentParticle().

◆ getMass()

Mdouble ParticleHandler::getMass ( ) const
606 {
607 #ifdef MERCURYDPM_USE_MPI
608  Mdouble massLocal = getMassLocal();
609  Mdouble massGlobal = 0.0;
610 
611  //Sum up over all domains
612  MPIContainer& communicator = MPIContainer::Instance();
613  communicator.allReduce(massLocal, massGlobal, MPI_SUM);
614 
615  return massGlobal;
616 #else
617  return getMassLocal();
618 #endif
619 }
Mdouble getMassLocal() const
Definition: ParticleHandler.cc:596

References getMassLocal(), and MPIContainer::Instance().

Referenced by BoundariesSelfTest::actionsAfterTimeStep(), FluxAndPeriodicBoundarySelfTest::actionsAfterTimeStep(), FluxBoundarySelfTest::actionsAfterTimeStep(), FixedClusterInsertionBoundary::checkBoundaryBeforeTimeStep(), RandomClusterInsertionBoundary::checkBoundaryBeforeTimeStep(), getCentreOfMass(), DPMBase::getTotalMass(), HeaterBoundaryTest::setupInitialConditions(), DPMBase::writeEneTimeStep(), and LawinenBox::writeEneTimeStep().

◆ getMassLocal()

Mdouble ParticleHandler::getMassLocal ( ) const
private
597 {
598  Mdouble m = 0;
599  for (auto p : *this)
600  if (!(p->isFixed() || p->isMPIParticle() || p->isPeriodicGhostParticle()))
601  m += p->getMass();
602  return m;
603 }

References m, and p.

Referenced by getMass().

◆ getMassTimesPosition()

Vec3D ParticleHandler::getMassTimesPosition ( ) const
631 {
632 #ifdef MERCURYDPM_USE_MPI
633  Vec3D massTimesPositionLocal = getMassTimesPositionLocal();
634  Vec3D massTimesPositionGlobal = {0.0, 0.0, 0.0};
635 
636  // Sum up over all domains
637  MPIContainer& communicator = MPIContainer::Instance();
638  communicator.allReduce(massTimesPositionLocal.X, massTimesPositionGlobal.X, MPI_SUM);
639  communicator.allReduce(massTimesPositionLocal.Y, massTimesPositionGlobal.Y, MPI_SUM);
640  communicator.allReduce(massTimesPositionLocal.Z, massTimesPositionGlobal.Z, MPI_SUM);
641 
642  return massTimesPositionGlobal;
643 #else
644  return getMassTimesPositionLocal();
645 #endif
646 }
Vec3D getMassTimesPositionLocal() const
Definition: ParticleHandler.cc:621
Mdouble Y
Definition: Kernel/Math/Vector.h:45
Mdouble Z
Definition: Kernel/Math/Vector.h:45
Mdouble X
the vector components
Definition: Kernel/Math/Vector.h:45

References getMassTimesPositionLocal(), MPIContainer::Instance(), Vec3D::X, Vec3D::Y, and Vec3D::Z.

Referenced by getCentreOfMass(), and DPMBase::writeEneTimeStep().

◆ getMassTimesPositionLocal()

Vec3D ParticleHandler::getMassTimesPositionLocal ( ) const
private
622 {
623  Vec3D com = {0, 0, 0};
624  for (auto p : *this)
625  if (!(p->isFixed() || p->isMPIParticle() || p->isPeriodicGhostParticle()))
626  com += p->getMass() * p->getPosition();
627  return com;
628 }

References p.

Referenced by getMassTimesPosition().

◆ getMeanRadius()

Mdouble ParticleHandler::getMeanRadius ( ) const
799 {
800 #ifdef MERCURYDPM_USE_MPI
801  Mdouble sumRadiusLocal = getSumRadiusLocal();
802  unsigned numberOfRealParticlesLocal = getNumberOfRealObjectsLocal();
803 
804  Mdouble sumRadiusGlobal = 0.0;
805  unsigned numberOfRealParticlesGlobal = 0;
806 
807  //Sum up over all domains
808  MPIContainer& communicator = MPIContainer::Instance();
809  communicator.allReduce(sumRadiusLocal, sumRadiusGlobal, MPI_SUM);
810  communicator.allReduce(numberOfRealParticlesLocal, numberOfRealParticlesGlobal, MPI_SUM);
811 
812  return sumRadiusGlobal/numberOfRealParticlesGlobal;
813 #else
814  return getSumRadiusLocal() / getSize();
815 #endif
816 }
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
Mdouble getSumRadiusLocal() const
Definition: ParticleHandler.cc:785

References getNumberOfRealObjectsLocal(), BaseHandler< BaseParticle >::getSize(), getSumRadiusLocal(), and MPIContainer::Instance().

Referenced by InitialConditions< SpeciesType >::getMeanRelativeContactRadius(), Sintering::getMeanRelativeContactRadius(), regimeForceUnitTest::getMeanRelativeContactRadius(), LaserOnLayer::LaserOnLayer(), GranuDrum::printTime(), Drum::printTime(), RotatingDrumWet::printTime(), RotatingDrumBidisperse::printTime(), regimeForceUnitTest::printTime(), and SingleParticle< SpeciesType >::writeEneTimeStep().

◆ getMomentum()

Vec3D ParticleHandler::getMomentum ( ) const
661 {
662  Vec3D momentum = {0, 0, 0};
663  for (auto p : *this)
664  if (!(p->isFixed() || p->isMPIParticle() || p->isPeriodicGhostParticle()))
665  momentum += p->getMomentum();
666  return getMPISum(momentum);
667 }

References getMPISum(), and p.

Referenced by DPMBase::getTotalMomentum(), ParticleBeam::printTime(), and LawinenBox::writeEneTimeStep().

◆ getName()

std::string ParticleHandler::getName ( ) const
overridevirtual

Returns the name of the handler, namely the string "ParticleHandler".

Returns
The string "ParticleHandler".

Implements BaseHandler< BaseParticle >.

1246 {
1247  return "ParticleHandler";
1248 }

◆ getNumberOfFixedObjects()

unsigned int ParticleHandler::getNumberOfFixedObjects ( ) const

Computes the number of fixed particles in the whole simulation.

1357 {
1358 #ifdef MERCURYDPM_USE_MPI
1359  unsigned int numberOfFixedParticlesLocal = getNumberOfFixedObjectsLocal();
1360  unsigned int numberOfFixedParticles = 0;
1361 
1362  MPIContainer& communicator = MPIContainer::Instance();
1363  communicator.allReduce(numberOfFixedParticlesLocal, numberOfFixedParticles, MPI_SUM);
1364  return numberOfFixedParticles;
1365 #else
1367 #endif
1368 }
unsigned int getNumberOfFixedObjectsLocal() const
Computes the number of Fixed particles on a local domain.
Definition: ParticleHandler.cc:1340

References getNumberOfFixedObjectsLocal(), and MPIContainer::Instance().

Referenced by LawinenBox::setupInitialConditions().

◆ getNumberOfFixedObjectsLocal()

unsigned int ParticleHandler::getNumberOfFixedObjectsLocal ( ) const

Computes the number of Fixed particles on a local domain.

Loops over all particles to check if the particle is fixed or not, in a local domain.

Returns
the number of fixed particles in a local domain
1341 {
1342  unsigned int numberOfFixedParticles = 0;
1343  for (BaseParticle* particle : *this)
1344  {
1345  if (particle->isFixed())
1346  {
1347  numberOfFixedParticles++;
1348  }
1349  }
1350  return numberOfFixedParticles;
1351 }
bool isFixed() const override
Is fixed Particle function. It returns whether a Particle is fixed or not, by checking its inverse Ma...
Definition: BaseParticle.h:72

References BaseParticle::isFixed().

Referenced by getNumberOfFixedObjects().

◆ getNumberOfFixedParticles()

unsigned int ParticleHandler::getNumberOfFixedParticles ( ) const

Gets the number of particles that are fixed.

1004 {
1005  return NFixedParticles_;
1006 }

References NFixedParticles_.

◆ getNumberOfObjects()

unsigned int ParticleHandler::getNumberOfObjects ( ) const
overridevirtual

Returns the number of objects in the container. In parallel code this practice is forbidden to avoid confusion with real and fake particles. If the size of the container is wished, call size()

Reimplemented from BaseHandler< BaseParticle >.

1324 {
1325 //#ifdef MERCURYDPM_USE_MPI
1326 // MPIContainer& communicator = MPIContainer::Instance();
1327 // if (communicator.getNumberOfProcessors() > 1)
1328 // {
1329 // logger(WARN,"When compiling with MPI please do not use getNumberOfObjects(). Instead use: getNumberOfRealObjectsLocal(), getNumberOfRealObjects() or getSize()");
1330 // }
1331 //#endif
1332  return getSize();
1333 }

References BaseHandler< BaseParticle >::getSize().

Referenced by Contact::actionsAfterSolve(), TwoByTwoMPIDomainMPI4Test::actionsAfterTimeStep(), GranularCollapse::actionsAfterTimeStep(), MercuryProblem::actionsAfterTimeStep(), LawinenBox::actionsBeforeTimeStep(), SmoothChute::actionsBeforeTimeStep(), AngleOfRepose::actionsBeforeTimeStep(), Chutebelt::actionsBeforeTimeStep(), ChutePeriodic::add_flow_particles(), SilbertPeriodic::add_flow_particles(), ChuteWithPeriodicInflow::AddContinuingBottom(), Chute::addFlowParticlesCompactly(), NautaMixer::addParticles(), statistics_while_running< T >::auto_set_domain(), statistics_while_running< T >::auto_set_z(), Mercury3Dclump::checkClumpForInteractionPeriodic(), DPMBase::checkParticleForInteractionLocalPeriodic(), ChuteWithContraction::ChuteWithContraction(), ChuteWithPeriodicInflowAndContinuingBottom::ChuteWithPeriodicInflowAndContinuingBottom(), ChuteWithPeriodicInflowAndContraction::ChuteWithPeriodicInflowAndContraction(), ChuteWithPeriodicInflowAndVariableBottom::ChuteWithPeriodicInflowAndVariableBottom(), ChuteWithPeriodicInflow::cleanChute(), ChuteWithContraction::cleanChute(), Funnel::cleanChute(), Chute::cleanChute(), ClosedCSCRestart::ClosedCSCRestart(), ClosedCSCRun::ClosedCSCRun(), ContractionWithPeriodicInflow::ContractionWithPeriodicInflow(), Slide::create_rough_wall(), AngleOfRepose::createBaseSpecies(), SilbertPeriodic::createBaseSpecies(), CSCInit::CSCInit(), ChuteWithPeriodicInflow::ExtendInWidth(), SphericalIndenter::getBedHeight(), getNumberOfUnfixedParticles(), getPSD(), InitialConditions< SpeciesType >::InitialConditions(), ChuteWithPeriodicInflow::integrateBeforeForceComputation(), ContactDetectionIntersectionOfWallsTest::introduceParticlesAtWall(), LawinenBox::LawinenBox(), main(), ChuteBottom::makeRoughBottom(), LiquidMigrationSelfTest::outputXBallsData(), SphericalIndenter::outputXBallsData(), LawinenBox::printTime(), GranuDrum::printTime(), GranuHeap::printTime(), ChuteWithPeriodicInflow::printTime(), Drum::printTime(), RotatingDrumWet::printTime(), VerticalMixer::printTime(), SilbertPeriodic::printTime(), vibratedBed::printTime(), RotatingDrumBidisperse::printTime(), Chute::printTime(), readAndCreateObject(), RotatingDrumBidisperse::RotatingDrumBidisperse(), CSCInit::save(), RotatingDrumBidisperseInitialise::save(), save(), ClosedCSCWalls::saveWalls(), CSCWalls::saveWalls(), MeshTriangle::setHandler(), LawinenBox::setupInitialConditions(), ClosedCSCWalls::setupInitialConditions(), CSCRun::setupInitialConditions(), CSCWalls::setupInitialConditions(), Binary::setupInitialConditions(), FreeCooling2DinWalls::setupInitialConditions(), FreeCooling2DinWallsDemo::setupInitialConditions(), FreeCooling3DDemoProblem::setupInitialConditions(), FreeCooling3DinWallsDemo::setupInitialConditions(), FreeCoolingDemoProblem::setupInitialConditions(), HourGlass2D::setupInitialConditions(), HourGlass::setupInitialConditions(), Cstatic2d::setupInitialConditions(), Cstatic3D::setupInitialConditions(), AngleOfRepose::setupInitialConditions(), SilbertPeriodic::setupInitialConditions(), InitialBed::setupInitialConditions(), Chutebelt::setupInitialConditions(), ParticleCreation::setupInitialConditions(), TriangulatedScrewSelfTest::setupInitialConditions(), TriangulatedStepWallSelfTest::setupInitialConditions(), TriangulatedWallSelfTest::setupInitialConditions(), SphericalIndenter::setupInitialConditions(), ScalingTestRun::setupInitialConditions(), GranularCollapse::setupInitialConditions(), Tutorial11::setupInitialConditions(), MD_demo::setupInitialConditions(), ChuteBottom::setupInitialConditions(), and SingleParticleIndenter::SingleParticleIndenter().

◆ getNumberOfRealObjects()

unsigned int ParticleHandler::getNumberOfRealObjects ( ) const

Returns the number of real objects (on all processors)

1303 {
1304 #ifdef MERCURYDPM_USE_MPI
1305  MPIContainer& communicator = MPIContainer::Instance();
1306  unsigned int numberOfRealParticles = getNumberOfRealObjectsLocal();
1307 
1308  // \todo MX: use an allreduce here
1309  //Combine the total number of Particles into one number on processor 0
1310  communicator.reduce(numberOfRealParticles, MPI_SUM);
1311 
1312  //Broadcast new number to all the processorsi
1313  communicator.broadcast(numberOfRealParticles);
1314  return numberOfRealParticles;
1315 #else
1316  return getSize();
1317 #endif
1318 }

References MPIContainer::broadcast(), getNumberOfRealObjectsLocal(), BaseHandler< BaseParticle >::getSize(), and MPIContainer::Instance().

Referenced by Membrane::createVertexParticles(), DPMBase::decompose(), getPSD(), and HourGlass2D::setupInitialConditions().

◆ getNumberOfRealObjectsLocal()

unsigned int ParticleHandler::getNumberOfRealObjectsLocal ( ) const

Returns the number of real objects on a local domain. MPI particles and periodic particles are neglected.

Returns
the number of real particles (neglecting MPI and periodic) on a local domain
Todo:
MX: in future also add the periodic mpi particles
1282 {
1283 #ifdef MERCURYDPM_USE_MPI
1284  const MPIContainer& communicator = MPIContainer::Instance();
1285  if (communicator.getNumberOfProcessors() > 1)
1286  {
1287  unsigned int numberOfFakeParticles = 0;
1288  numberOfFakeParticles += getDPMBase()->domainHandler.getCurrentDomain()->getNumberOfTrueMPIParticles();
1290  logger.assert_debug(numberOfFakeParticles <= getSize(), "More fake particles than getSize()");
1291  return (getSize() - numberOfFakeParticles);
1292  }
1293  else
1294  {
1295  return getSize();
1296  }
1297 #else
1298  return getSize();
1299 #endif
1300 }
PeriodicBoundaryHandler periodicBoundaryHandler
Internal handler that deals with periodic boundaries, especially in a parallel build.
Definition: DPMBase.h:1463
DomainHandler domainHandler
An object of the class DomainHandler which deals with parallel code.
Definition: DPMBase.h:1468
Domain * getCurrentDomain()
Gets the domain assigned to the processor.
Definition: DomainHandler.cc:228
unsigned int getNumberOfTrueMPIParticles()
Obtains the number of particles in the particleHandler that are MPIParticles, but NOT periodic partic...
Definition: Domain.cc:1658
std::size_t getNumberOfProcessors() const
Get the total number of processors participating in this simulation.
Definition: MpiContainer.cc:83
unsigned int getNumberOfPeriodicGhostParticles()
Returns the number of particles that are flagged is periodicGhostParticle.
Definition: PeriodicBoundaryHandler.cc:335

References DPMBase::domainHandler, DomainHandler::getCurrentDomain(), BaseHandler< BaseParticle >::getDPMBase(), PeriodicBoundaryHandler::getNumberOfPeriodicGhostParticles(), MPIContainer::getNumberOfProcessors(), Domain::getNumberOfTrueMPIParticles(), BaseHandler< BaseParticle >::getSize(), MPIContainer::Instance(), logger, and DPMBase::periodicBoundaryHandler.

Referenced by MaserRepeatedOutInMPI2Test::actionsAfterSolve(), TwoByTwoMPIDomainMPI4Test::actionsAfterTimeStep(), getMeanRadius(), getNumberOfRealObjects(), and DPMBase::outputXBallsData().

◆ getNumberOfUnfixedParticles()

unsigned int ParticleHandler::getNumberOfUnfixedParticles ( ) const

Gets the number of particles that are not fixed.

1012 {
1014 }
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

References getNumberOfObjects(), and NFixedParticles_.

◆ getParticleAttribute()

template<typename T >
std::enable_if<std::is_scalar<T>::value, T>::type ParticleHandler::getParticleAttribute ( std::function< T(BaseParticle *)>  attribute,
AttributeType  type 
) const
inline

Computes an attribute type (min/max/..) of a particle attribute(position/velocity) in the global domain.

Many functions the particleHandler return a pointer to a BaseParticle, i.e. the fastest particle, however in parallel this is not useful as pointers can't be send across processes. This function gives the flexibility to find a global particle extremum such as the fastest particle Velocity or lowest particle position.

Parameters
[in]attributeA function that obtains a scalar from a particle. i.e. position or radius
[in]typeAn AttributeType tells this function what to do with the attribute, take the maximum or the minimum
Returns
Returns a scalar value representing the AttributeType of the Attribute: In easier terms returns the max/min/... of a particle attribute such as radius
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  }
@ MIN
Definition: ParticleHandler.h:18
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
default
Definition: calibrate.py:45

References calibrate::default, ERROR, getParticleAttributeLocal(), MPIContainer::Instance(), logger, MAX, MIN, and compute_granudrum_aor::type.

◆ getParticleAttributeLocal()

template<typename T >
std::enable_if<std::is_scalar<T>::value, T>::type ParticleHandler::getParticleAttributeLocal ( std::function< T(BaseParticle *)>  attribute,
AttributeType  type 
) const
inline

Computes an attribute type (min/max/..) of a particle attribute (position/velocity) in a local domain.

Computes an attribute type (min/max/..) of a particle attribute (position/velocity) in a local domain

Many functions the particleHandler return a pointer to a BaseParticle, i.e. the fastest particle, however in parallel this is not useful as pointers can't be send across processes. This function gives the flexibility to find a global particle extremum such as the fastest particle Velocity or lowest particle position.

Parameters
[in]attributeA function that obtains a scalar from a particle. i.e. position or radius
[in]typeAn AttributeType tells this function what to do with the attribute, take the maximum or the minimum
Returns
Returns a scalar value representing the AttributeType of the Attribute: In easier terms returns the max/min/... of a particle attribute such as radius
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  }

References calibrate::default, ERROR, logger, MAX, MIN, BaseHandler< BaseParticle >::objects_, compute_granudrum_aor::type, and WARN.

Referenced by getParticleAttribute().

◆ getPSD()

PSD ParticleHandler::getPSD ( std::vector< Mdouble bins = {},
Mdouble  scaleFactor = 1.0 
) const

Gets the PSD from the particles currently in the handler.

When no bins are provided, they will be created linearly between the smallest and largest particle radius in the handler. When bins are provided and don't include the smallest/largest particle size currently in the handler, they will be added. For the largest particle size this is a necessity, otherwise there is no way to account for larger particles. For the smallest particle size this isn't strictly necessary, however it is preferred to have the PSD start with zero probability.

Parameters
binsOptionally, the bins of interest.
scaleFactorOptionally, value to scale the bins by, to convert them to radii. E.g. 0.5 when diameters are provided.
Returns
1428 {
1429  // check that there is at least one particle in the handler
1430  if (getNumberOfObjects()==0) {
1431  logger(WARN, "getPSD: Empty PSD is returned, as particleHandler is empty.");
1432  return PSD();
1433  }
1434 
1435  // check that the bins vector is set
1436 #ifdef MERCURYDPM_USE_MPI
1437  double rMin = getSmallestInteractionRadius();
1438  double rMax = getLargestInteractionRadius();
1439 #else
1440  double rMin = getSmallestParticle()->getRadius();
1441  double rMax = getLargestParticle()->getRadius();
1442 #endif
1443  if (bins.empty()) {
1444  // Number of bins needed with average of 20 particles in each bin.
1445  int n = (int)std::ceil(getNumberOfObjects()/20);
1446  // Limit number of bins to a certain maximum and minimum.
1447  n = std::min(100, n);
1448  n = std::max(5, n);
1449  bins.reserve(n);
1450  double dr = (rMax-rMin)/(n-1);
1451  for (int i = 0; i<n-1; i++) {
1452  bins.push_back(rMin + i * dr);
1453  }
1454  } else {
1455  // make sure values are in increasing order
1456  std::sort(bins.begin(), bins.end());
1457  // and convert to radius
1458  for (auto& b : bins) {
1459  b *= scaleFactor;
1460  }
1461  }
1462 
1463  // set the first/last bin to include the smallest/largest particle (first bin to make PSD start with zero probability)
1464  if (bins.front() > rMin) {
1465  bins.insert(bins.begin(), rMin);
1466  }
1467  if (bins.back() < rMax) {
1468  bins.push_back(rMax);
1469  }
1470 
1471  // PSD vector, all probabilities initialised to 0
1472  std::vector<DistributionElements> des;
1473  des.reserve(bins.size());
1474  for (const Mdouble bin : bins) {
1475  des.push_back({ bin, 0.0 });
1476  }
1477 
1478  // calculate the probabilities
1479  for (const auto& p : *this) {
1480  if (p->isFixed() || p->isMPIParticle() || p->isPeriodicGhostParticle()) {
1481  continue;
1482  }
1483  // the first "bin"/radius is skipped, since the probabilities for each bin are assigned to its upper radius
1484  // i.e. the PSD should start with zero probability
1485  for (auto bin = des.begin() + 1; bin != des.end(); bin++) {
1486  if (p->getRadius() <= bin->internalVariable) {
1487  bin->probability += 1.0 / getNumberOfRealObjects();
1488  break;
1489  }
1490  }
1491  }
1492 
1493  // create and return PSD
1494  PSD psd;
1496  return psd;
1497 }
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Scalar * b
Definition: benchVecAdd.cpp:17
Mdouble getRadius() const
Returns the particle's radius.
Definition: BaseParticle.h:331
Contains a vector with radii and probabilities of a user defined particle size distribution (PSD)
Definition: PSD.h:47
@ PROBABILITYDENSITY_NUMBER_DISTRIBUTION
void setPSDFromVector(std::vector< DistributionElements > psd, TYPE PSDType)
Sets the PSD from a vector with DistributionElements.
Definition: PSD.cc:513
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 getNumberOfRealObjects() const
Returns the number of real objects (on all processors)
Definition: ParticleHandler.cc:1302
BaseParticle * getSmallestParticle() const
Gets a pointer to the smallest BaseParticle (by interactionRadius) in this ParticleHandler of the loc...
Definition: ParticleHandler.cc:505
return int(ret)+1
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 ceil(const bfloat16 &a)
Definition: BFloat16.h:644
bins
Definition: UniformPSDSelfTest.py:19

References b, UniformPSDSelfTest::bins, Eigen::bfloat16_impl::ceil(), getLargestInteractionRadius(), getLargestParticle(), getNumberOfObjects(), getNumberOfRealObjects(), BaseParticle::getRadius(), getSmallestInteractionRadius(), getSmallestParticle(), i, int(), logger, max, min, n, p, PSD::PROBABILITYDENSITY_NUMBER_DISTRIBUTION, PSD::setPSDFromVector(), and WARN.

Referenced by saveNumberPSDtoCSV().

◆ getRotationalEnergy()

Mdouble ParticleHandler::getRotationalEnergy ( ) const
581 {
582 #ifdef MERCURYDPM_USE_MPI
583  Mdouble rotationalEnergyLocal = getRotationalEnergyLocal();
584  Mdouble rotationalEnergyGlobal = 0.0;
585 
586  //sum up over all domains
587  MPIContainer& communicator = MPIContainer::Instance();
588  communicator.allReduce(rotationalEnergyLocal, rotationalEnergyGlobal, MPI_SUM);
589 
590  return rotationalEnergyGlobal;
591 #else
592  return getRotationalEnergyLocal();
593 #endif
594 }
Mdouble getRotationalEnergyLocal() const
Definition: ParticleHandler.cc:567

References getRotationalEnergyLocal(), and MPIContainer::Instance().

Referenced by DPMBase::writeEneTimeStep(), and LawinenBox::writeEneTimeStep().

◆ getRotationalEnergyLocal()

Mdouble ParticleHandler::getRotationalEnergyLocal ( ) const
private
568 {
569  Mdouble ene = 0;
570  for (auto p : *this)
571  {
572  if (!(p->isMPIParticle() || p->isPeriodicGhostParticle()))
573  {
574  ene += p->getRotationalEnergy();
575  }
576  }
577  return ene;
578 }

References p.

Referenced by getRotationalEnergy().

◆ getSmallestInteractionRadius()

Mdouble ParticleHandler::getSmallestInteractionRadius ( ) const

Returns the smallest interaction radius.

731 {
732 #ifdef MERCURYDPM_USE_MPI
733  //Compute the local value
734  Mdouble smallestInteractionRadiusLocal = getSmallestInteractionRadiusLocal();
735  Mdouble smallestInteractionRadiusGlobal = 0.0;
736 
737  //Obtain the global value
738  MPIContainer& communicator = MPIContainer::Instance();
739  communicator.allReduce(smallestInteractionRadiusLocal, smallestInteractionRadiusGlobal, MPI_MIN);
740 
741  return smallestInteractionRadiusGlobal;
742 
743 #else
745 #endif
746 }
Mdouble getSmallestInteractionRadiusLocal() const
Returns the smallest interaction radius of the current domain.
Definition: ParticleHandler.cc:715

References getSmallestInteractionRadiusLocal(), and MPIContainer::Instance().

Referenced by Contact::actionsAfterSolve(), getPSD(), and main().

◆ getSmallestInteractionRadiusLocal()

Mdouble ParticleHandler::getSmallestInteractionRadiusLocal ( ) const

Returns the smallest interaction radius of the current domain.

716 {
717  if (!(getSmallestParticleLocal() == nullptr))
718  {
720  }
721  else
722  {
723  return 0.0;
724  }
725 }
BaseParticle * getSmallestParticleLocal() const
Gets a pointer to the smallest BaseParticle (by interactionRadius) in this ParticleHandler of the loc...
Definition: ParticleHandler.cc:496

References BaseParticle::getMaxInteractionRadius(), and getSmallestParticleLocal().

Referenced by MercuryBase::getHGridTargetMinInteractionRadius(), and getSmallestInteractionRadius().

◆ getSmallestParticle()

BaseParticle * ParticleHandler::getSmallestParticle ( ) const

Gets a pointer to the smallest BaseParticle (by interactionRadius) in this ParticleHandler of the local domain. When mercury is running in parallel this function throws an error since a pointer to another domain is useless.

Returns
A pointer to the to the smallest BaseParticle (by interactionRadius) in this ParticleHandler.
506 {
507 #ifdef MERCURYDPM_USE_MPI
508  logger(WARN,"getSmallestParticle should not be used in parallel; use getSmallestInteractionRadius or "
509  "ParticleSpecies::getSmallestParticleMass instead");
510  return nullptr;
511 #else
512  return getSmallestParticleLocal();
513 #endif
514 }

References getSmallestParticleLocal(), logger, and WARN.

Referenced by RotatingDrumBidisperseInitialise::actionsAfterTimeStep(), getPSD(), and CFDDEMCoupleTest::setupInitialConditions().

◆ getSmallestParticleLocal()

BaseParticle * ParticleHandler::getSmallestParticleLocal ( ) const

Gets a pointer to the smallest BaseParticle (by interactionRadius) in this ParticleHandler of the local domain.

Returns
A pointer to the to the smallest BaseParticle (by interactionRadius) in this ParticleHandler.
497 {
498  return smallestParticle_;
499 }

References smallestParticle_.

Referenced by getSmallestInteractionRadiusLocal(), getSmallestParticle(), RotatingDrum::setupInitialConditions(), and DPMBase::writeFstatHeader().

◆ getSumRadiusLocal()

Mdouble ParticleHandler::getSumRadiusLocal ( ) const
private
Returns
The sum of all particle radii
786 {
787  Mdouble sumRadius = 0;
788  for (BaseParticle* const p : objects_)
789  {
790  if (!(p->isMPIParticle() || p->isPeriodicGhostParticle()))
791  {
792  sumRadius += p->getRadius();
793  }
794  }
795  return sumRadius;
796 }

References BaseHandler< BaseParticle >::objects_, and p.

Referenced by getMeanRadius().

◆ getVolume()

Mdouble ParticleHandler::getVolume ( ) const
1262 {
1263 #ifdef MERCURYDPM_USE_MPI
1264  Mdouble volumeLocal = getVolumeLocal();
1265  Mdouble volumeGlobal = 0.0;
1266 
1267  // sum up over all domains
1268  MPIContainer& communicator = MPIContainer::Instance();
1269  communicator.allReduce(volumeLocal, volumeGlobal, MPI_SUM);
1270 
1271  return volumeGlobal;
1272 #else
1273  return getVolumeLocal();
1274 #endif
1275 }
Mdouble getVolumeLocal() const
Definition: ParticleHandler.cc:1250

References getVolumeLocal(), and MPIContainer::Instance().

Referenced by PSDManualInsertionSelfTest::actionsAfterSolve(), ShearStage::actionsAfterTimeStep(), LeesEdwardsDemo::actionsAfterTimeStep(), ShiftingConstantMassFlowMaserBoundarySelfTest::actionsAfterTimeStep(), ShiftingMaserBoundarySelfTest::actionsAfterTimeStep(), TimeDependentPeriodicBoundary3DSelfTest::actionsAfterTimeStep(), TimeDependentPeriodicBoundaryTest::actionsAfterTimeStep(), FixedClusterInsertionBoundary::checkBoundaryBeforeTimeStep(), RandomClusterInsertionBoundary::checkBoundaryBeforeTimeStep(), and GranuHeap::GranuHeap().

◆ getVolumeLocal()

Mdouble ParticleHandler::getVolumeLocal ( ) const
private
1251 {
1252  Mdouble volume = 0;
1253  for (auto p : *this)
1254  {
1255  if (!(p->isFixed() || p->isPeriodicGhostParticle() || p->isMPIParticle()))
1256  volume += p->getVolume();
1257  }
1258  return volume;
1259 }

References p.

Referenced by getVolume().

◆ operator=()

ParticleHandler & ParticleHandler::operator= ( const ParticleHandler rhs)

Assignment operator.

Parameters
[in]rhsThe ParticleHandler on the right hand side of the assignment.

This is not a copy assignment operator! It copies the DPMBase and all BaseParticle, and sets the other variables to 0. After that, it computes the smallest and largest particle in this handler.

65 {
66  if (this != &rhs)
67  {
68  clear();
69  largestParticle_ = nullptr;
70  smallestParticle_ = nullptr;
72  if (!objects_.empty())
73  {
76  }
77  }
78  logger(DEBUG, "ParticleHandler::operator = (const ParticleHandler& rhs) finished");
79  return *this;
80 }

References clear(), computeLargestParticle(), computeSmallestParticle(), BaseHandler< BaseParticle >::copyContentsFromOtherHandler(), DEBUG, largestParticle_, logger, BaseHandler< BaseParticle >::objects_, and smallestParticle_.

◆ readAndAddObject()

void ParticleHandler::readAndAddObject ( std::istream &  is)
overridevirtual
     \brief Create a new particle, based on the information from old-style restart data.
    &zwj;/

void readAndCreateOldObject(std::istream &is, const std::string& type);

/*!
   \brief Create a new particle in the WallHandler, based on the information provided in a restart file.
Parameters
[in]isThe input stream from which the information is read.

Implements BaseHandler< BaseParticle >.

1093 {
1095  addExistingObject(o);
1096 }
void addExistingObject(BaseParticle *P) override
Adds a BaseParticle to the ParticleHandler.
Definition: ParticleHandler.cc:110
BaseParticle * readAndCreateObject(std::istream &is)
Create a new particle, based on the information provided in a restart file.
Definition: ParticleHandler.cc:1063

References addExistingObject(), and readAndCreateObject().

Referenced by DPMBase::read().

◆ readAndCreateObject()

BaseParticle * ParticleHandler::readAndCreateObject ( std::istream &  is)

Create a new particle, based on the information provided in a restart file.

Parameters
[in]isThe input stream from which the information is read.
Todo:
check what is special about this case
1064 {
1065  std::string type;
1066  BaseParticle* particle = nullptr;
1067  if (getNumberOfObjects()>0) logger(DEBUG,"size % % %",getNumberOfObjects(),getLastObject()->getVelocity(),getLastObject()->getRadius());
1068  if (getStorageCapacity() > 0)
1069  {
1070  is >> type;
1071  logger(DEBUG, "ParticleHandler::readAndCreateObject(is): type %", type);
1072  particle = createObject(type);
1073  particle->setHandler(this);
1074  particle->read(is);
1075  }
1076  else //for insertion boundaries
1077  {
1079  is >> type;
1080  logger(DEBUG, "ParticleHandler::readAndCreateObject(is): type %", type);
1081  particle = createObject(type);
1082  particle->setHandler(this);
1083  particle->read(is);
1084  }
1085  particle->setSpecies(getDPMBase()->speciesHandler.getObject(particle->getIndSpecies()));
1086  return particle;
1087 }
unsigned int getStorageCapacity() const
Gets the storage capacity of this BaseHandler.
Definition: BaseHandler.h:670
BaseParticle * getLastObject()
Gets a pointer to the last Object in this BaseHandler.
Definition: BaseHandler.h:642
unsigned int getIndSpecies() const
Returns the index of the species associated with the interactable object.
Definition: BaseInteractable.h:67
void setHandler(ParticleHandler *handler)
Sets the pointer to the particle's ParticleHandler.
Definition: BaseParticle.cc:640
void read(std::istream &is) override
Particle read function, which accepts an std::istream as input.
Definition: BaseParticle.cc:368
virtual void setSpecies(const ParticleSpecies *species)
Definition: BaseParticle.cc:798
static BaseParticle * createObject(const std::string &type)
Reads BaseParticle into the ParticleHandler from restart data.
Definition: ParticleHandler.cc:1019
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286

References createObject(), DEBUG, BaseHandler< BaseParticle >::getDPMBase(), BaseInteractable::getIndSpecies(), BaseHandler< BaseParticle >::getLastObject(), getNumberOfObjects(), BaseHandler< BaseParticle >::getStorageCapacity(), logger, BaseParticle::read(), BaseParticle::setHandler(), BaseParticle::setSpecies(), oomph::Global_string_for_annotation::string(), and compute_granudrum_aor::type.

Referenced by InsertionBoundary::read(), and readAndAddObject().

◆ removedFixedParticle()

void ParticleHandler::removedFixedParticle ( )

Decrement of the number of fixed particles.

Todo:
MX: For Jonny, is this still required, keeping the parallel code in mind?
1238 {
1239  NFixedParticles_--;
1240 }

References NFixedParticles_.

Referenced by BaseParticle::unfix(), and BaseParticle::~BaseParticle().

◆ removeGhostObject()

void ParticleHandler::removeGhostObject ( unsigned int  index)

Removes a BaseParticle from the ParticleHandler without a global check, this is only to be done for mpi routines.

Parameters
[in]indexThe index of which BaseParticle has to be removed from this ParticleHandler.

The BaseParticle with index is removed and the last BaseParticle in the vector is moved to its position. It also removes the BaseParticle from the HGrid. This function is only called by the parallel routines where we are sure the particles have been flushed out of the communication boundaries.

417 {
418 #ifdef MERCURYDPM_USE_MPI
419 #ifdef CONTACT_LIST_HGRID
420  getDPMBase()->getPossibleContactList().remove_ParticlePosibleContacts(getObject(id));
421 #endif
423  // broadcast the particle removal
424  getDPMBase()->handleParticleRemoval(getObject(index)->getId());
426 #else
427  logger(ERROR, "This function should only be used interally for mpi routines");
428 #endif
429 }
virtual void removeObject(unsigned const int index)
Removes an Object from the BaseHandler.
Definition: BaseHandler.h:453
BaseParticle * getObject(const unsigned int id)
Gets a pointer to the Object at the specified index in the BaseHandler.
Definition: BaseHandler.h:621
virtual void hGridRemoveParticle(BaseParticle *obj UNUSED)
Definition: DPMBase.cc:1700

References ERROR, BaseHandler< BaseParticle >::getDPMBase(), BaseHandler< BaseParticle >::getObject(), DPMBase::handleParticleRemoval(), DPMBase::hGridRemoveParticle(), logger, and BaseHandler< T >::removeObject().

Referenced by DeletionBoundary::checkBoundaryAfterParticleMoved(), DPMBase::checkInteractionWithBoundaries(), PeriodicBoundaryHandler::clearCommunicationLists(), and DPMBase::deleteGhostParticles().

◆ removeLastObject()

void ParticleHandler::removeLastObject ( )

Removes the last BaseParticle from the ParticleHandler.

Function that removes the last object from this ParticleHandler. It also removes the particle from the HGrid.

436 {
437 #ifdef CONTACT_LIST_HGRID
438  getDPMBase()->getPossibleContactList().remove_ParticlePosibleContacts(getLastObject());
439 #endif
443 }
void removeLastObject()
Removes the last Object from the BaseHandler.
Definition: BaseHandler.h:519

References BaseHandler< BaseParticle >::getDPMBase(), BaseHandler< BaseParticle >::getLastObject(), DPMBase::handleParticleRemoval(), DPMBase::hGridRemoveParticle(), and BaseHandler< T >::removeLastObject().

Referenced by Penetration::setupInitialConditions().

◆ removeObject()

void ParticleHandler::removeObject ( unsigned int  index)
overridevirtual

Removes a BaseParticle from the ParticleHandler.

Parameters
[in]indexThe index of which BaseParticle has to be removed from this ParticleHandler.

The BaseParticle with index is removed and the last BaseParticle in the vector is moved to its position. It also removes the BaseParticle from the HGrid. When running this in parallel a warning is thrown that deleting particles might cause havoc in the communication between boundaries, as the deleted particle might still be in one of the boundary lists

Reimplemented from BaseHandler< BaseParticle >.

389 {
390 #ifdef MERCURYDPM_USE_MPI
391  MPIContainer& communicator = MPIContainer::Instance();
392  if (communicator.getNumberOfProcessors() > 1 )
393  {
394  logger(WARN, "[ParticleHandler::removeObject(const unsigned int)] Using the function removeObject in parallel could lead to out-of-sync communication, Instead use deletion boundaries");
395  }
396 #endif
397 #ifdef CONTACT_LIST_HGRID
398  getDPMBase()->getPossibleContactList().remove_ParticlePosibleContacts(getObject(id));
399 #endif
401  // broadcast the particle removal
402  getDPMBase()->handleParticleRemoval(getObject(index)->getId());
404 }

References BaseHandler< BaseParticle >::getDPMBase(), MPIContainer::getNumberOfProcessors(), BaseHandler< BaseParticle >::getObject(), DPMBase::handleParticleRemoval(), DPMBase::hGridRemoveParticle(), MPIContainer::Instance(), logger, BaseHandler< T >::removeObject(), and WARN.

Referenced by AngleOfRepose::actionsBeforeTimeStep(), PeriodicWallsWithSlidingFrictionUnitTest::actionsBeforeTimeStep(), SilbertPeriodic::add_flow_particles(), CircularPeriodicBoundary::checkBoundaryAfterParticleMoved(), ConstantMassFlowMaserBoundary::checkBoundaryAfterParticleMoved(), DeletionBoundary::checkBoundaryAfterParticleMoved(), ChuteWithContraction::ChuteWithContraction(), ChuteWithPeriodicInflowAndContraction::ChuteWithPeriodicInflowAndContraction(), ChuteWithPeriodicInflow::cleanChute(), ChuteWithContraction::cleanChute(), Funnel::cleanChute(), Chute::cleanChute(), ContactDetectionTester::cleanup(), ContractionWithPeriodicInflow::ContractionWithPeriodicInflow(), ChuteBottom::makeRoughBottom(), CurvyChute::recreateBottom(), DPMBase::removeDuplicatePeriodicParticles(), CSCWalls::saveWalls(), and AngleOfRepose::setupInitialConditions().

◆ saveNumberPSDtoCSV()

void ParticleHandler::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_NUMBER_DISTRIBUTION.

Parameters
csvFileNameThe filename to write to.
diameterBinsOptionally, the bins of interest.
1401 {
1402  // get the PSD, with scale factor 0.5 to scale diameter to radius
1403  PSD psd = getPSD(diameterBins, 0.5);
1404  // get the PSD vector by the required type, with scale factor 2 to scale radius to diameter
1405  std::vector<DistributionElements> psdVector = psd.getParticleSizeDistributionByType(PSD::TYPE::PROBABILITYDENSITY_NUMBER_DISTRIBUTION, 2.0);
1406 
1407  logger(INFO,"Writing PSD to %",csvFileName);
1408  std::ofstream csv(csvFileName);
1409  logger.assert_always(csv.is_open(),"File % could not be opened",csvFileName);
1410  csv << "Diameter" << ',' << "Probability Number Distribution" << '\n';
1411  for (const auto& bin : psdVector) {
1412  csv << bin.internalVariable << ',' << bin.probability << '\n';
1413  }
1414  csv.close();
1415 }
std::vector< DistributionElements > getParticleSizeDistributionByType(TYPE psdType, Mdouble scalingFactor=1.0) const
Gets the PSD vector, converted to the preferred type.
Definition: PSD.cc:1067
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

References PSD::getParticleSizeDistributionByType(), getPSD(), INFO, logger, and PSD::PROBABILITYDENSITY_NUMBER_DISTRIBUTION.

Referenced by main().

◆ write()

Member Data Documentation

◆ largestParticle_

BaseParticle* ParticleHandler::largestParticle_
private

A pointer to the largest BaseParticle (by interactionRadius) in this ParticleHandler.

Todo:
TW: note that checkExtrema gets called if a particle gets created and its Species and Radius gets set, even if it's not yet included in the particleHandler! This is necessary to check a not included particle for overlaps before inserting it into the handler. Not sure if this is a sensible structure; to be discussed. Note: these statistic now include mpi and ghost particles as well

Referenced by checkExtrema(), checkExtremaOnDelete(), clear(), computeLargestParticle(), getLargestParticleLocal(), operator=(), ParticleHandler(), and ~ParticleHandler().

◆ NFixedParticles_

unsigned int ParticleHandler::NFixedParticles_
private

◆ smallestParticle_

BaseParticle* ParticleHandler::smallestParticle_
private

A pointer to the smallest BaseParticle (by interactionRadius) in this ParticleHandler.

Referenced by checkExtrema(), checkExtremaOnDelete(), clear(), computeSmallestParticle(), getSmallestParticleLocal(), operator=(), ParticleHandler(), and ~ParticleHandler().


The documentation for this class was generated from the following files: