BaseParticle Class Referenceabstract

#include <BaseParticle.h>

+ Inheritance diagram for BaseParticle:

Public Member Functions

 BaseParticle ()
 Basic Particle constructor, creates an Particle at (0,0,0) with radius, mass and inertia equal to 1. More...
 
 BaseParticle (const BaseParticle &p)
 Particle copy constructor, which accepts as input a reference to a Particle. It creates a copy of this Particle and all it's information. Usually it is better to use the copy() function for polymorphism. More...
 
 BaseParticle (const ParticleSpecies *s)
 
 ~BaseParticle () override
 Particle destructor, needs to be implemented and checked if it removes tangential spring information. More...
 
virtual BaseParticlecopy () const =0
 Particle copy method. It calls to copy constructor of this Particle, useful for polymorphism. More...
 
virtual Mdouble getVolume () const
 Get Particle volume function, which required a reference to the Species vector. It returns the volume of the Particle. More...
 
void fixParticle ()
 Fix Particle function. It fixes a Particle by setting its inverse mass and inertia and velocities to zero. More...
 
bool isFixed () const override
 Is fixed Particle function. It returns whether a Particle is fixed or not, by checking its inverse Mass. More...
 
bool isMPIParticle () const
 Indicates if this particle is a ghost in the MPI domain. More...
 
void setMPIParticle (bool flag)
 Flags the mpi particle status. More...
 
bool isInMPIDomain ()
 Indicates if the particle is in the communication zone of the mpi domain. More...
 
void setInMPIDomain (bool flag)
 Flags the status of the particle if wether it is in the communication zone or not. More...
 
bool isInPeriodicDomain () const
 Indicates if the particle is in the periodic boundary communication zone. More...
 
void setInPeriodicDomain (bool flag)
 Flags the status of the particle whether it is in the periodic communication zone or not. More...
 
bool isPeriodicGhostParticle () const
 Indicates if this particle is a ghost in the periodic boundary. More...
 
void setPeriodicGhostParticle (bool flag)
 Flags the status of the particle to be a ghost in periodic boundary or not. More...
 
bool isMaserParticle () const
 Indicates if this particle belongs to the maser boundary. More...
 
void setMaserParticle (bool flag)
 Flags the status of the particle if it belongs to the maser boundary or not. More...
 
void setCommunicationComplexity (unsigned complexity)
 Set the communication complexity of the particle. More...
 
unsigned getCommunicationComplexity ()
 Obtains the communication complexity of the particle. More...
 
void setPeriodicComplexity (std::vector< int > complexity)
 Set the periodic communication complexity of the particle. More...
 
void setPeriodicComplexity (int index, int value)
 Set the periodic communication complexity of the particle. More...
 
const std::vector< int > & getPeriodicComplexity ()
 Obtains the periodic communication complexity of the particle. More...
 
void setPreviousPeriodicComplexity (std::vector< int > complexity)
 Set the previous periodic communication complexity of the paritcle. More...
 
const std::vector< int > & getPreviousPeriodicComplexity () const
 Sets the previous periodic communication complexity of the particle. More...
 
int getPeriodicComplexity (int index)
 Gets the periodic communication complexity of a certain boundary. More...
 
void unfix ()
 Unfix Particle function, which required a reference to the Species vector. It unfixes a Particle by computing the Particles mass and inertia. More...
 
void read (std::istream &is) override
 Particle read function, which accepts an std::istream as input. More...
 
virtual void oldRead (std::istream &is)
 
void write (std::ostream &os) const override
 Particle print function, which accepts an std::ostream as input. More...
 
std::string getName () const override
 Returns the name of the object. More...
 
virtual void setInfo (Mdouble info)
 Sets some user-defined information about this object (by default, species ID). More...
 
virtual Mdouble getInfo () const
 Returns some user-defined information about this object (by default, species ID). More...
 
void setTimeStamp (unsigned timeStamp)
 
unsigned getTimeStamp () const
 
void printHGrid (std::ostream &os) const
 Adds particle's HGrid level and cell coordinates to an ostream. More...
 
unsigned int getHGridLevel () const
 Returns particle's HGrid level. More...
 
BaseParticlegetHGridNextObject () const
 Returns pointer to next object in particle's HGrid level & cell. More...
 
BaseParticlegetHGridPrevObject () const
 Returns pointer to previous object in particle's HGrid level & cell. More...
 
int getHGridX () const
 Returns particle's HGrid cell X-coordinate. More...
 
int getHGridY () const
 Returns particle's HGrid cell Y-coordinate. More...
 
int getHGridZ () const
 Returns particle's HGrid cell Z-coordinate. More...
 
MatrixSymmetric3D getInvInertia () const
 Returns the inverse of the particle's inertia tensor. More...
 
Mdouble getInvMass () const override
 Returns the inverse of the particle's mass. More...
 
Mdouble getCurvature (const Vec3D &labFixedCoordinates) const override
 
virtual Mdouble getKineticEnergy () const
 Calculates the particle's translational kinetic energy. More...
 
virtual Mdouble getRotationalEnergy () const
 Calculates the particle's rotational kinetic energy. More...
 
Mdouble getGravitationalEnergy () const
 Calculates the particle's gravitational energy. More...
 
Mdouble getMass () const
 Returns the particle's mass. More...
 
virtual Mdouble getSurfaceArea () const
 
Vec3D getMomentum () const
 
MatrixSymmetric3D getInertia () const
 
Vec3D getAngularMomentum () const
 
BaseParticlegetPeriodicFromParticle () const
 Returns the 'original' particle this one's a periodic copy of. More...
 
Mdouble getRadius () const
 Returns the particle's radius. More...
 
Mdouble getMaxInteractionRadius () const
 Returns the particle's interaction radius, which might be different from radius_ (e.g., when dealing with wet particles) More...
 
Mdouble getInteractionDistance (const BaseInteractable *i) const
 Returns the interactionDistance_ of the mixed species of this particle and the particle or wall i. More...
 
Mdouble getSumOfInteractionRadii (const BaseParticle *particle) const
 returns the sum of the radii plus the interactionDistance More...
 
Mdouble getWallInteractionRadius (const BaseWall *wall) const
 returns the radius plus the interactionDistance More...
 
const Vec3DgetPreviousPosition () const
 Returns the particle's position in the previous time step. More...
 
const Vec3D getDisplacement2 (Mdouble xmin, Mdouble xmax, Mdouble ymin, Mdouble ymax, Mdouble zmin, Mdouble zmax, Mdouble t) const
 
virtual void setInertia ()
 
void setInertia (MatrixSymmetric3D inertia)
 Sets the particle's inertia_ (and adjusts invInertia_ accordingly) More...
 
void setInverseInertia (MatrixSymmetric3D inverseInertia)
 Sets the particle's inertia_ (and adjusts invInertia_ accordingly) More...
 
void setInfiniteInertia ()
 Sets the particle's inertia_ to 'infinite' (1e20) and its invInertia_ to 0. More...
 
void setPeriodicFromParticle (BaseParticle *p)
 Assigns the pointer to the 'original' particle this one's a periodic copy of (used in periodic boundary condition implementations). More...
 
void setHGridX (const int x)
 Sets the particle's HGrid cell X-coordinate. More...
 
void setHGridY (const int y)
 Sets the particle's HGrid cell Y-coordinate. More...
 
void setHGridZ (const int z)
 Sets the particle's HGrid cell Z-coordinate. More...
 
void setHGridLevel (const unsigned int level)
 Sets the particle's HGrid level. More...
 
void setHGridNextObject (BaseParticle *p)
 Sets the pointer to the next object in the particle's HGrid cell & level. More...
 
void setHGridPrevObject (BaseParticle *p)
 Sets the pointer to the previous object in the particle's HGrid cell & level. More...
 
virtual void setRadius (Mdouble radius)
 Sets the particle's radius_ (and adjusts the mass_ accordingly, based on the particle's species) More...
 
virtual Vec3D getAxes () const
 Only ustilised in case of superquadric particles. Had to create a virtual function to allow function access in writeVTK function in the particle handler. More...
 
virtual Mdouble getExponentEps1 () const
 Only ustilised in case of superquadric particles. Had to create a virtual function to allow function access in writeVTK function in the particle handler. More...
 
virtual Mdouble getExponentEps2 () const
 Only ustilised in case of superquadric particles. Had to create a virtual function to allow function access in writeVTK function in the particle handler. More...
 
virtual void setAxes (const Vec3D &axes)
 Only ustilised in case of superquadric particles. More...
 
virtual void setExponents (const Mdouble &eps1, const Mdouble &eps2)
 Only ustilised in case of superquadric particles. More...
 
MERCURYDPM_DEPRECATED void setMass (Mdouble mass)
 Sets the particle's mass. More...
 
void setMassForP3Statistics (Mdouble mass)
 Sets the particle's mass This function should not be used, but is necessary to extend the CG toolbox to non-spherical particles. More...
 
void setPreviousPosition (const Vec3D &pos)
 Sets the particle's position in the previous time step. More...
 
void movePrevious (const Vec3D &posMove)
 Adds a vector to the particle's previousPosition_. More...
 
void accelerate (const Vec3D &vel)
 Increases the particle's velocity_ by the given vector. More...
 
void angularAccelerate (const Vec3D &angVel)
 Increases the particle's angularVelocity_ by the given vector. More...
 
void setHandler (ParticleHandler *handler)
 Sets the pointer to the particle's ParticleHandler. More...
 
ParticleHandlergetHandler () const
 Returns pointer to the particle's ParticleHandler. More...
 
BaseInteractiongetInteractionWith (BaseParticle *P, unsigned timeStamp, InteractionHandler *interactionHandler) override
 Checks if particle is in interaction with given particle P, and if so, returns vector of pointer to the associated BaseInteraction object (else returns empty vector). More...
 
virtual bool isInContactWith (const BaseParticle *P) const
 Get whether or not this particle is in contact with the given particle. More...
 
virtual void integrateBeforeForceComputation (double time, double timeStep)
 First step of Velocity Verlet integration. More...
 
virtual void integrateAfterForceComputation (double time, double timeStep)
 Second step of Velocity Verlet integration. More...
 
unsigned int getParticleDimensions () const
 Returns the particle's dimensions (either 2 or 3). More...
 
MERCURYDPM_DEPRECATED void setIndSpecies (unsigned int indSpecies) override
 
virtual void setSpecies (const ParticleSpecies *species)
 
virtual void actionsBeforeTimeStep ()
 
virtual unsigned getNumberOfFieldsVTK () const
 
virtual std::string getTypeVTK (unsigned i) const
 
virtual std::string getNameVTK (unsigned i) const
 
virtual std::vector< MdoublegetFieldVTK (unsigned i) const
 
virtual void actionsAfterTimeStep ()
 
virtual bool isSphericalParticle () const =0
 
const HGridCellgetHGridCell () const
 
virtual void computeMass (const ParticleSpecies &s)
 Computes the particle's (inverse) mass and inertia. More...
 
BaseParticlegetClump () const
 
bool isClump () const
 Checks if particle is a clump (container) More...
 
bool isPebble () const
 Checks if particle is a pebble (belongs to a clump) More...
 
virtual Vec3D getCenterOfMass ()
 
virtual void actionsAfterAddObject ()
 
- Public Member Functions inherited from BaseInteractable
 BaseInteractable ()
 Default BaseInteractable constructor. More...
 
 BaseInteractable (const BaseInteractable &p)
 Copy constructor. More...
 
 ~BaseInteractable () override
 Destructor, it simply destructs the BaseInteractable and all the objects it contains. More...
 
unsigned int getIndSpecies () const
 Returns the index of the species associated with the interactable object. More...
 
const ParticleSpeciesgetSpecies () const
 Returns a pointer to the species of this BaseInteractable. More...
 
void setSpecies (const ParticleSpecies *species)
 Sets the species of this BaseInteractable. More...
 
const Vec3DgetForce () const
 Returns the force on this BaseInteractable. More...
 
const Vec3DgetTorque () const
 Returns the torque on this BaseInteractable. More...
 
void setForce (const Vec3D &force)
 Sets the force on this BaseInteractable. More...
 
void setTorque (const Vec3D &torque)
 Sets the torque on this BaseInteractable. More...
 
void addForce (const Vec3D &addForce)
 Adds an amount to the force on this BaseInteractable. More...
 
void addTorque (const Vec3D &addTorque)
 Adds an amount to the torque on this BaseInteractable. More...
 
virtual void resetForceTorque (int numberOfOMPthreads)
 
void sumForceTorqueOMP ()
 
const Vec3DgetPosition () const
 Returns the position of this BaseInteractable. More...
 
const QuaterniongetOrientation () const
 Returns the orientation of this BaseInteractable. More...
 
virtual void setPosition (const Vec3D &position)
 Sets the position of this BaseInteractable. More...
 
void setOrientationViaNormal (Vec3D normal)
 Sets the orientation of this BaseInteractable by defining the vector that results from the rotation of the (1,0,0) vector. More...
 
void setOrientationViaEuler (Vec3D eulerAngle)
 Sets the orientation of this BaseInteractable by defining the euler angles. More...
 
virtual void setOrientation (const Quaternion &orientation)
 Sets the orientation of this BaseInteractable. More...
 
virtual void move (const Vec3D &move)
 Moves this BaseInteractable by adding an amount to the position. More...
 
virtual void rotate (const Vec3D &angularVelocityDt)
 Rotates this BaseInteractable. More...
 
const std::vector< BaseInteraction * > & getInteractions () const
 Returns a list of interactions which belong to this interactable. More...
 
void addInteraction (BaseInteraction *I)
 Adds an interaction to this BaseInteractable. More...
 
bool removeInteraction (BaseInteraction *I)
 Removes an interaction from this BaseInteractable. More...
 
void copyInteractionsForPeriodicParticles (const BaseInteractable &p)
 Copies interactions to this BaseInteractable whenever a periodic copy made. More...
 
void setVelocity (const Vec3D &velocity)
 set the velocity of the BaseInteractable. More...
 
void setAngularVelocity (const Vec3D &angularVelocity)
 set the angular velocity of the BaseInteractble. More...
 
void addVelocity (const Vec3D &velocity)
 adds an increment to the velocity. More...
 
void addAngularVelocity (const Vec3D &angularVelocity)
 add an increment to the angular velocity. More...
 
virtual const Vec3DgetVelocity () const
 Returns the velocity of this interactable. More...
 
virtual const Vec3DgetAngularVelocity () const
 Returns the angular velocity of this interactable. More...
 
void setPrescribedPosition (const std::function< Vec3D(double)> &prescribedPosition)
 Allows the position of an infinite mass interactable to be prescribed. More...
 
void applyPrescribedPosition (double time)
 Computes the position from the user defined prescribed position function. More...
 
void setPrescribedVelocity (const std::function< Vec3D(double)> &prescribedVelocity)
 Allows the velocity of an infinite mass interactable to be prescribed. More...
 
void applyPrescribedVelocity (double time)
 Computes the velocity from the user defined prescribed velocity function. More...
 
void setPrescribedOrientation (const std::function< Quaternion(double)> &prescribedOrientation)
 Allows the orientation of the infinite mass interactbale to be prescribed. More...
 
void applyPrescribedOrientation (double time)
 Computes the orientation from the user defined prescribed orientation function. More...
 
void setPrescribedAngularVelocity (const std::function< Vec3D(double)> &prescribedAngularVelocity)
 Allows the angular velocity of the infinite mass interactable to be prescribed. More...
 
void applyPrescribedAngularVelocity (double time)
 Computes the angular velocity from the user defined prescribed angular velocity. More...
 
virtual const Vec3D getVelocityAtContact (const Vec3D &contact) const
 Returns the velocity at the contact point, use by many force laws. More...
 
void integrateBeforeForceComputation (double time, double timeStep)
 This is part of integrate routine for objects with infinite mass. More...
 
void integrateAfterForceComputation (double time, double timeStep)
 This is part of the integration routine for objects with infinite mass. More...
 
virtual bool isFaceContact (const Vec3D &normal) const
 
- Public Member Functions inherited from BaseObject
 BaseObject ()=default
 Default constructor. More...
 
 BaseObject (const BaseObject &p)=default
 Copy constructor, copies all the objects BaseObject contains. More...
 
virtual ~BaseObject ()=default
 virtual destructor More...
 
virtual void moveInHandler (unsigned int index)
 Except that it is virtual, it does the same thing as setIndex() does. More...
 
void setIndex (unsigned int index)
 Allows one to assign an index to an object in the handler/container. More...
 
void setId (unsigned long id)
 Assigns a unique identifier to each object in the handler (container) which remains constant even after the object is deleted from the container/handler. More...
 
unsigned int getIndex () const
 Returns the index of the object in the handler. More...
 
unsigned int getId () const
 Returns the unique identifier of any particular object. More...
 
void setGroupId (unsigned groupId)
 
unsigned getGroupId () const
 

Public Attributes

Mdouble radius_
 
Mdouble invMass_
 Particle radius_. More...
 
MatrixSymmetric3D invInertia_
 Inverse Particle mass (for computation optimization) More...
 
BaseParticleclumpParticle_
 Function that updates necessary quantities of a clump particle after adding a pebble. More...
 
bool isPebble_
 pointer to a clump particle (for a pebble) More...
 
bool isClump_
 The particle is pebble. More...
 

Private Attributes

ParticleHandlerhandler_
 Inverse Particle inverse inertia (for computation optimization) More...
 
HGridCell hGridCell
 
BaseParticlehGridNextObject_
 
BaseParticlehGridPrevObject_
 Pointer to the next Particle in the same HGrid cell. More...
 
BaseParticleperiodicFromParticle_
 Pointer to the previous Particle in the same HGrid cell. More...
 
bool isMPIParticle_
 Pointer to originating Particle. More...
 
bool isInMPIDomain_
 returns true if the particle acts as an MPI particle instead of a real particle More...
 
unsigned communicationComplexity_
 returns true if it flagged as being in MPI domain More...
 
bool isInPeriodicDomain_
 
bool isPeriodicGhostParticle_
 bool that indicates if a particle is in the periodic domain of any boundary More...
 
std::vector< intpreviousPeriodicComplexity_
 Indicates if the particle is a ghost particle of a periodic particle. More...
 
std::vector< intperiodicComplexity_
 Indicates the periodic complexity at previous time step. More...
 
bool isMaserParticle_
 Indicates the periodic complexity at current time step. Used to update periodic status. More...
 
Vec3D previousPosition_
 Indicates if this particle belongs to the maser boundary or is released into the wide open world. More...
 
Mdouble info_
 
unsigned timeStamp_
 

Friends

void ParticleSpecies::computeMass (BaseParticle *) const
 Particle's position at previous time step. More...
 

Detailed Description

Since r3648, BaseParticle is an abstract class. Use SphericalParticle for a 'basic' particle.

Constructor & Destructor Documentation

◆ BaseParticle() [1/3]

BaseParticle::BaseParticle ( )

Basic Particle constructor, creates an Particle at (0,0,0) with radius, mass and inertia equal to 1.

default constructor, creates an Particle at (0,0,0) with radius, mass and inertia equal to 1

13 {
14  handler_ = nullptr;
15  radius_ = 1.0;
16  invMass_ = 1.0;
17  invInertia_ = MatrixSymmetric3D(1, 0, 0, 1, 0, 1);
18 
20 
21  periodicFromParticle_ = nullptr;
22  isMPIParticle_ = false;
23  isInMPIDomain_ = false;
24  isInPeriodicDomain_ = false;
26  isMaserParticle_ = false;
28  periodicComplexity_ = std::vector<int>(0);
29  previousPeriodicComplexity_ = std::vector<int>(0);
30 #ifdef CONTACT_LIST_HGRID
31  firstPossibleContact = nullptr;
32 #endif
33  hGridNextObject_ = nullptr;
34  hGridPrevObject_ = nullptr;
35 
36  hGridCell.setHGridLevel(99999);
37  hGridCell.setHGridX(99999);
38  hGridCell.setHGridY(99999);
39  hGridCell.setHGridZ(99999);
40 
41  info_ = std::numeric_limits<double>::quiet_NaN();
42  timeStamp_ = 0;
43  isClump_ = false;
44  isPebble_ = false;
45  logger(DEBUG, "BaseParticle::BaseParticle() finished");
46 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ DEBUG
BaseInteractable()
Default BaseInteractable constructor.
Definition: BaseInteractable.cc:20
const Vec3D & getPosition() const
Returns the position of this BaseInteractable.
Definition: BaseInteractable.h:197
bool isMaserParticle_
Indicates the periodic complexity at current time step. Used to update periodic status.
Definition: BaseParticle.h:688
bool isClump_
The particle is pebble.
Definition: BaseParticle.h:717
BaseParticle * periodicFromParticle_
Pointer to the previous Particle in the same HGrid cell.
Definition: BaseParticle.h:676
unsigned timeStamp_
Definition: BaseParticle.h:703
Vec3D previousPosition_
Indicates if this particle belongs to the maser boundary or is released into the wide open world.
Definition: BaseParticle.h:690
bool isInPeriodicDomain_
Definition: BaseParticle.h:684
BaseParticle * hGridPrevObject_
Pointer to the next Particle in the same HGrid cell.
Definition: BaseParticle.h:673
std::vector< int > periodicComplexity_
Indicates the periodic complexity at previous time step.
Definition: BaseParticle.h:687
unsigned communicationComplexity_
returns true if it flagged as being in MPI domain
Definition: BaseParticle.h:681
Mdouble radius_
Definition: BaseParticle.h:650
std::vector< int > previousPeriodicComplexity_
Indicates if the particle is a ghost particle of a periodic particle.
Definition: BaseParticle.h:686
bool isPebble_
pointer to a clump particle (for a pebble)
Definition: BaseParticle.h:715
MatrixSymmetric3D invInertia_
Inverse Particle mass (for computation optimization)
Definition: BaseParticle.h:652
Mdouble invMass_
Particle radius_.
Definition: BaseParticle.h:651
bool isMPIParticle_
Pointer to originating Particle.
Definition: BaseParticle.h:679
Mdouble info_
Definition: BaseParticle.h:698
BaseParticle * hGridNextObject_
Definition: BaseParticle.h:672
HGridCell hGridCell
Definition: BaseParticle.h:670
bool isPeriodicGhostParticle_
bool that indicates if a particle is in the periodic domain of any boundary
Definition: BaseParticle.h:685
ParticleHandler * handler_
Inverse Particle inverse inertia (for computation optimization)
Definition: BaseParticle.h:659
bool isInMPIDomain_
returns true if the particle acts as an MPI particle instead of a real particle
Definition: BaseParticle.h:680
void setHGridZ(int HGridZ)
Definition: HGridCell.h:60
void setHGridX(int HGridX)
Definition: HGridCell.h:40
void setHGridY(int HGridY)
Definition: HGridCell.h:50
void setHGridLevel(unsigned int HGridLevel)
Definition: HGridCell.h:70
Implementation of a 3D symmetric matrix.
Definition: MatrixSymmetric.h:16

References communicationComplexity_, DEBUG, BaseInteractable::getPosition(), handler_, hGridCell, hGridNextObject_, hGridPrevObject_, info_, invInertia_, invMass_, isClump_, isInMPIDomain_, isInPeriodicDomain_, isMaserParticle_, isMPIParticle_, isPebble_, isPeriodicGhostParticle_, logger, periodicComplexity_, periodicFromParticle_, previousPeriodicComplexity_, previousPosition_, radius_, HGridCell::setHGridLevel(), HGridCell::setHGridX(), HGridCell::setHGridY(), HGridCell::setHGridZ(), and timeStamp_.

◆ BaseParticle() [2/3]

BaseParticle::BaseParticle ( const BaseParticle p)

Particle copy constructor, which accepts as input a reference to a Particle. It creates a copy of this Particle and all it's information. Usually it is better to use the copy() function for polymorphism.

Constructor that copies most of the properties of the given particle. Please note that not everything is copied, for example the position in the HGrid is not determined yet by the end of this constructor. It also does not copy the interactions and the pointer to the handler that handles this particle. Use with care.

Parameters
[in,out]pReference to the BaseParticle this one should become a copy of.
58 {
59  handler_ = nullptr;
60  radius_ = p.radius_;
61  invMass_ = p.getInvMass();
62  invInertia_ = p.getInvInertia();
63 
64  previousPosition_ = p.previousPosition_;
65 
66  hGridNextObject_ = nullptr;
67  hGridPrevObject_ = nullptr;
68 
69  hGridCell.setHGridLevel(p.getHGridLevel());
70  hGridCell.setHGridX(99999);
71  hGridCell.setHGridY(99999);
72  hGridCell.setHGridZ(99999);
73 
74  periodicFromParticle_ = p.periodicFromParticle_;
75  isMPIParticle_ = p.isMPIParticle_;
76  isInMPIDomain_ = p.isInMPIDomain_;
77  isInPeriodicDomain_ = p.isInPeriodicDomain_;
78  isMaserParticle_ = p.isMaserParticle_;
79  isPeriodicGhostParticle_ = p.isPeriodicGhostParticle_;
80  communicationComplexity_ = p.communicationComplexity_;
81  isClump_ = p.isClump();
82  isPebble_ = p.isPebble();
83  //periodicComplexity_ = p.periodicComplexity_;
84  //previousPeriodicComplexity_ = p.previousPeriodicComplexity_;
85 #ifdef CONTACT_LIST_HGRID
86  firstPossibleContact = nullptr;
87 #endif
88 
89  info_ = p.info_;
90  timeStamp_ = p.timeStamp_;
91  logger(DEBUG, "BaseParticle::BaseParticle(BaseParticle &p) finished");
92 }
float * p
Definition: Tutorial_Map_using.cpp:9

References communicationComplexity_, DEBUG, handler_, hGridCell, hGridNextObject_, hGridPrevObject_, info_, invInertia_, invMass_, isClump_, isInMPIDomain_, isInPeriodicDomain_, isMaserParticle_, isMPIParticle_, isPebble_, isPeriodicGhostParticle_, logger, p, periodicFromParticle_, previousPosition_, radius_, HGridCell::setHGridLevel(), HGridCell::setHGridX(), HGridCell::setHGridY(), HGridCell::setHGridZ(), and timeStamp_.

◆ BaseParticle() [3/3]

BaseParticle::BaseParticle ( const ParticleSpecies s)
explicit
95  : BaseParticle()
96 {
97  setSpecies(s);
98 #ifdef CONTACT_LIST_HGRID
99  firstPossibleContact = nullptr;
100 #endif
101  logger(DEBUG, "BaseParticle::BaseParticle(BaseSpecies &s) finished");
102 }
virtual void setSpecies(const ParticleSpecies *species)
Definition: BaseParticle.cc:798
BaseParticle()
Basic Particle constructor, creates an Particle at (0,0,0) with radius, mass and inertia equal to 1.
Definition: BaseParticle.cc:11
RealScalar s
Definition: level1_cplx_impl.h:130

References DEBUG, logger, s, and setSpecies().

◆ ~BaseParticle()

BaseParticle::~BaseParticle ( )
override

Particle destructor, needs to be implemented and checked if it removes tangential spring information.

Destructor. It asks the ParticleHandler to check if this was the smallest or largest particle and adjust itself accordingly.

109 {
110 
111  if (getHandler() != nullptr)
112  {
114  if (isFixed())
116  }
117  logger(DEBUG, "BaseParticle::~BaseParticle() of particle % finished.", getId());
118 
119 }
unsigned int getId() const
Returns the unique identifier of any particular object.
Definition: BaseObject.h:104
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
ParticleHandler * getHandler() const
Returns pointer to the particle's ParticleHandler.
Definition: BaseParticle.cc:650
void checkExtremaOnDelete(BaseParticle *P)
Checks if the extrema of this ParticleHandler needs updating when a particle is deleted.
Definition: ParticleHandler.cc:1197
void removedFixedParticle()
Decrement of the number of fixed particles.
Definition: ParticleHandler.cc:1237

References ParticleHandler::checkExtremaOnDelete(), DEBUG, getHandler(), BaseObject::getId(), isFixed(), logger, and ParticleHandler::removedFixedParticle().

Member Function Documentation

◆ accelerate()

void BaseParticle::accelerate ( const Vec3D vel)

Increases the particle's velocity_ by the given vector.

increases the the particle's velocity_ (BaseInteractable member) by adding the given vector.

Parameters
[in]velvector to be added to the velocity_
621 {
622  addVelocity(vel);
623 }
void addVelocity(const Vec3D &velocity)
adds an increment to the velocity.
Definition: BaseInteractable.h:291

References BaseInteractable::addVelocity().

Referenced by integrateAfterForceComputation(), ClumpParticle::integrateAfterForceComputation(), integrateBeforeForceComputation(), ClumpParticle::integrateBeforeForceComputation(), and MovingIntersectionOfWallsUnitTest_MovingReferenceFrame::setupInitialConditions().

◆ actionsAfterAddObject()

virtual void BaseParticle::actionsAfterAddObject ( )
inlinevirtual

Methods and attributes necessary for clumped particles

Reimplemented in ClumpParticle.

711 {}

◆ actionsAfterTimeStep()

virtual void BaseParticle::actionsAfterTimeStep ( )
inlinevirtual
608 {}

◆ actionsBeforeTimeStep()

virtual void BaseParticle::actionsBeforeTimeStep ( )
inlinevirtual
598 {}

◆ angularAccelerate()

void BaseParticle::angularAccelerate ( const Vec3D angVel)

Increases the particle's angularVelocity_ by the given vector.

increases the particle's angularVelocity_ (BaseInteractable member) by adding the given vector.

Parameters
[in]angVelvector to be added to the angularVelocity_
631 {
632  addAngularVelocity(angVel);
633 }
void addAngularVelocity(const Vec3D &angularVelocity)
add an increment to the angular velocity.
Definition: BaseInteractable.cc:348

References BaseInteractable::addAngularVelocity().

Referenced by integrateAfterForceComputation(), and integrateBeforeForceComputation().

◆ computeMass()

void BaseParticle::computeMass ( const ParticleSpecies s)
virtual

Computes the particle's (inverse) mass and inertia.

Reimplemented in SuperQuadricParticle, and ClumpParticle.

856  {
857  if (isFixed()) return;
858  if (getParticleDimensions()==3) {
859  invMass_ = 1.0 / (4.0 / 3.0 * constants::pi * getRadius() * getRadius() * getRadius() * s.getDensity());
860  invInertia_ = MatrixSymmetric3D(1, 0, 0, 1, 0, 1) / (.4 * getMass() * mathsFunc::square(getRadius()));
861  } else {
862  invMass_ = 1.0 / (constants::pi * getRadius() * getRadius() * s.getDensity());
863  invInertia_ = MatrixSymmetric3D(1, 0, 0, 1, 0, 1) / (.5 * getMass() * mathsFunc::square(getRadius()));
864  }
865 };
unsigned int getParticleDimensions() const
Returns the particle's dimensions (either 2 or 3).
Definition: BaseParticle.cc:764
Mdouble getRadius() const
Returns the particle's radius.
Definition: BaseParticle.h:331
Mdouble getMass() const
Returns the particle's mass.
Definition: BaseParticle.h:305
const Mdouble pi
Definition: ExtendedMath.h:23
T square(const T val)
squares a number
Definition: ExtendedMath.h:86

References getMass(), getParticleDimensions(), getRadius(), invInertia_, invMass_, isFixed(), constants::pi, s, and mathsFunc::square().

Referenced by ChuteWithContraction::create_inflow_particle().

◆ copy()

◆ fixParticle()

void BaseParticle::fixParticle ( )

Fix Particle function. It fixes a Particle by setting its inverse mass and inertia and velocities to zero.

Fixes a BaseParticle by setting its inverse mass and inertia and velocities to zero.

151 {
152 // //
153 // Clump f;
154 // f.invInertia2_ = MatrixSymmetric3D(0, 0, 0, 0, 0, 0);
155  //
156  invMass_ = 0.0;
157  invInertia_ = MatrixSymmetric3D(0, 0, 0, 0, 0, 0);
158  setVelocity(Vec3D(0.0, 0.0, 0.0));
159  setAngularVelocity(Vec3D(0.0, 0.0, 0.0));
160  if (getHandler())
162 }
void setAngularVelocity(const Vec3D &angularVelocity)
set the angular velocity of the BaseInteractble.
Definition: BaseInteractable.cc:338
void setVelocity(const Vec3D &velocity)
set the velocity of the BaseInteractable.
Definition: BaseInteractable.cc:328
void addedFixedParticle()
Increment of the number of fixed particles.
Definition: ParticleHandler.cc:1232
Definition: Kernel/Math/Vector.h:30

References ParticleHandler::addedFixedParticle(), getHandler(), invInertia_, invMass_, BaseInteractable::setAngularVelocity(), and BaseInteractable::setVelocity().

Referenced by LawinenBox::LawinenBox(), ParticleBeam::ParticleBeam(), DPMBase::setFixedParticles(), Cstatic2d::setupInitialConditions(), Cstatic3D::setupInitialConditions(), and SphericalIndenter::setupInitialConditions().

◆ getAngularMomentum()

Vec3D BaseParticle::getAngularMomentum ( ) const
575 {
577 }
virtual const Vec3D & getAngularVelocity() const
Returns the angular velocity of this interactable.
Definition: BaseInteractable.cc:319
virtual const Vec3D & getVelocity() const
Returns the velocity of this interactable.
Definition: BaseInteractable.cc:307
static MatrixSymmetric3D inverse(const MatrixSymmetric3D &A)
Computes the inverse of a matrix; exits if the inverse doesn't exist.
Definition: MatrixSymmetric.cc:267
static Vec3D cross(const Vec3D &a, const Vec3D &b)
Calculates the cross product of two Vec3D: .
Definition: Vector.cc:143

References Vec3D::cross(), BaseInteractable::getAngularVelocity(), BaseInteractable::getPosition(), BaseInteractable::getVelocity(), MatrixSymmetric3D::inverse(), invInertia_, and invMass_.

◆ getAxes()

Vec3D BaseParticle::getAxes ( ) const
virtual

Only ustilised in case of superquadric particles. Had to create a virtual function to allow function access in writeVTK function in the particle handler.

Reimplemented in SuperQuadricParticle.

837 { return Vec3D(0, 0, 0); }

◆ getCenterOfMass()

virtual Vec3D BaseParticle::getCenterOfMass ( )
inlinevirtual
644 {return Vec3D(0,0,0);}

◆ getClump()

BaseParticle* BaseParticle::getClump ( ) const
inline
623  {
624  return clumpParticle_;
625  }
BaseParticle * clumpParticle_
Function that updates necessary quantities of a clump particle after adding a pebble.
Definition: BaseParticle.h:713

References clumpParticle_.

Referenced by Mercury3Dclump::computeInternalForce().

◆ getCommunicationComplexity()

unsigned BaseParticle::getCommunicationComplexity ( )

Obtains the communication complexity of the particle.

185 {
187 }

References communicationComplexity_.

Referenced by Domain::updateParticles().

◆ getCurvature()

Mdouble BaseParticle::getCurvature ( const Vec3D labFixedCoordinates) const
inlineoverridevirtual

returns the inverse radius, or curvature, of the surface. This value is zero for walls and gets overridden for particles that have finite radius

Todo:
should be wall-type dependent

Reimplemented from BaseInteractable.

Reimplemented in SuperQuadricParticle.

280  { return 1.0/radius_; }

References radius_.

◆ getDisplacement2()

const Vec3D BaseParticle::getDisplacement2 ( Mdouble  xmin,
Mdouble  xmax,
Mdouble  ymin,
Mdouble  ymax,
Mdouble  zmin,
Mdouble  zmax,
Mdouble  t 
) const
Todo:
see .cc file. \TWH
Todo:
Rewrite, redefine (TW). Is only used in StatisticsVector.hcc, consider moving to that class.
470 {
472  if (xmax > xmin && fabs(disp.X) > .5 * (xmax - xmin))
473  {
474  if (disp.X > 0)
475  disp.X -= xmax - xmin;
476  else
477  disp.X += xmax - xmin;
478  }
479  if (ymax > ymin && fabs(disp.Y) > .5 * (ymax - ymin))
480  {
481  if (disp.Y > 0)
482  disp.Y -= ymax - ymin;
483  else
484  disp.Y += ymax - ymin;
485  }
486  if (zmax > zmin && fabs(disp.Z) > .5 * (zmax - zmin))
487  {
488  if (disp.Z > 0)
489  disp.Z -= zmax - zmin;
490  else
491  disp.Z += zmax - zmin;
492  }
493  disp /= t;
494  return disp;
495 }
const Vec3D & getPreviousPosition() const
Returns the particle's position in the previous time step.
Definition: BaseParticle.h:378
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
double zmax
Definition: time_harmonic_fourier_decomposed_linear_elasticity/cylinder/cylinder.cc:71
double zmin
Definition: time_harmonic_fourier_decomposed_linear_elasticity/cylinder/cylinder.cc:69
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117
t
Definition: plotPSD.py:36

References boost::multiprecision::fabs(), BaseInteractable::getPosition(), getPreviousPosition(), plotPSD::t, Vec3D::X, Vec3D::Y, Vec3D::Z, Global_Parameters::zmax, and Global_Parameters::zmin.

◆ getExponentEps1()

double BaseParticle::getExponentEps1 ( ) const
virtual

Only ustilised in case of superquadric particles. Had to create a virtual function to allow function access in writeVTK function in the particle handler.

Reimplemented in SuperQuadricParticle.

840 { return 0; }

◆ getExponentEps2()

double BaseParticle::getExponentEps2 ( ) const
virtual

Only ustilised in case of superquadric particles. Had to create a virtual function to allow function access in writeVTK function in the particle handler.

Reimplemented in SuperQuadricParticle.

843 { return 0; }

◆ getFieldVTK()

std::vector< Mdouble > BaseParticle::getFieldVTK ( unsigned  i) const
virtual

Reimplemented in ClumpParticle.

832 {
833  return std::vector<Mdouble>();
834 }

◆ getGravitationalEnergy()

Mdouble BaseParticle::getGravitationalEnergy ( ) const

Calculates the particle's gravitational energy.

Gravitational energy is the potential energy stored in teh particles position due to the gravity field. This is a relative term, so we need to define what zero energy means: The gravitational energy of a particle is zero when its center of mass is at the origin.

Returns
the particle's gravitational energy

◆ getHandler()

◆ getHGridCell()

◆ getHGridLevel()

◆ getHGridNextObject()

BaseParticle* BaseParticle::getHGridNextObject ( ) const
inline

Returns pointer to next object in particle's HGrid level & cell.

Returns the next object in the particle's HGrid cell

Returns
pointer to the next object in the particle's HGrid cell
226  { return hGridNextObject_; }

References hGridNextObject_.

Referenced by Mercury3D::hGridFindContactsWithinTargetCell(), Mercury2D::hGridFindContactsWithinTargetCell(), Mercury2D::hGridRemoveParticle(), and Mercury3D::hGridRemoveParticle().

◆ getHGridPrevObject()

BaseParticle* BaseParticle::getHGridPrevObject ( ) const
inline

Returns pointer to previous object in particle's HGrid level & cell.

Returns the previous object in the particle's HGrid cell

Returns
pointer to the previous object in the particle's HGrid cell
234  { return hGridPrevObject_; }

References hGridPrevObject_.

Referenced by Mercury2D::hGridRemoveParticle(), and Mercury3D::hGridRemoveParticle().

◆ getHGridX()

int BaseParticle::getHGridX ( ) const
inline

Returns particle's HGrid cell X-coordinate.

Returns
the particle's HGrid cell's X-coordinate
249  { return hGridCell.getHGridX(); }
int getHGridX() const
Definition: HGridCell.h:35

References HGridCell::getHGridX(), and hGridCell.

Referenced by Mercury2D::computeInternalForces(), Mercury3D::computeInternalForces(), Mercury2D::hGridRemoveParticle(), Mercury2D::hGridUpdateParticle(), and Mercury3D::hGridUpdateParticle().

◆ getHGridY()

int BaseParticle::getHGridY ( ) const
inline

Returns particle's HGrid cell Y-coordinate.

Returns
the particle's HGrid cell's Y-coordinate
256  { return hGridCell.getHGridY(); }
int getHGridY() const
Definition: HGridCell.h:45

References HGridCell::getHGridY(), and hGridCell.

Referenced by Mercury2D::computeInternalForces(), Mercury3D::computeInternalForces(), Mercury2D::hGridRemoveParticle(), Mercury2D::hGridUpdateParticle(), and Mercury3D::hGridUpdateParticle().

◆ getHGridZ()

int BaseParticle::getHGridZ ( ) const
inline

Returns particle's HGrid cell Z-coordinate.

Returns
the particle's HGrid cell's Z-coordinate
263  { return hGridCell.getHGridZ(); }
int getHGridZ() const
Definition: HGridCell.h:55

References HGridCell::getHGridZ(), and hGridCell.

Referenced by Mercury3D::computeInternalForces(), and Mercury3D::hGridUpdateParticle().

◆ getInertia()

◆ getInfo()

Mdouble BaseParticle::getInfo ( ) const
virtual

Returns some user-defined information about this object (by default, species ID).

343 {
344  if (std::isnan(info_))
345  return getSpecies()->getId();
346  else
347  return info_;
348 }
const ParticleSpecies * getSpecies() const
Returns a pointer to the species of this BaseInteractable.
Definition: BaseInteractable.h:87
#define isnan(X)
Definition: main.h:109

References BaseObject::getId(), BaseInteractable::getSpecies(), info_, and isnan.

◆ getInteractionDistance()

Mdouble BaseParticle::getInteractionDistance ( const BaseInteractable i) const
inline

Returns the interactionDistance_ of the mixed species of this particle and the particle or wall i.

353  {
354  //const auto mixedSpecies = getSpecies()->getHandler()->getMixedObject(getSpecies(),particle->getSpecies());
355  //return mixedSpecies->getInteractionDistance();
356  return getSpecies()->getMixedSpecies(i->getSpecies())->getInteractionDistance();
357  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Mdouble getInteractionDistance() const
returns the largest separation distance at which adhesive short-range forces can occur.
Definition: BaseSpecies.h:125
const BaseSpecies * getMixedSpecies(const ParticleSpecies *s) const
Definition: ParticleSpecies.cc:218

References BaseSpecies::getInteractionDistance(), ParticleSpecies::getMixedSpecies(), BaseInteractable::getSpecies(), and i.

Referenced by getSumOfInteractionRadii(), and getWallInteractionRadius().

◆ getInteractionWith()

BaseInteraction * BaseParticle::getInteractionWith ( BaseParticle P,
unsigned  timeStamp,
InteractionHandler interactionHandler 
)
overridevirtual

Checks if particle is in interaction with given particle P, and if so, returns vector of pointer to the associated BaseInteraction object (else returns empty vector).

Creates/updates a BaseInteraction object, treating the interaction between this particle and a given one, in case there is an overlap between the two.

Parameters
[in]Pparticle to check the interaction with
[in]timeStamptime stamp to be assigned to the interaction object (i.e., the current time)
[in,out]interactionHandlerBaseInteraction container from where the interaction is retrieved, and to which it is assigned (if it is a new interaction).
Returns
the pointer to the interaction object (if the particles overlap), or 0 (if they don't overlap).

Implements BaseInteractable.

Reimplemented in SuperQuadricParticle.

669 {
670  //get the normal (from P away from the contact)
671  const Vec3D branchVector = P->getPosition() - getPosition();
672  //Get the square of the distance between particle i and particle j
673  const Mdouble distanceSquared = Vec3D::getLengthSquared(branchVector);
674  //const auto species = interactionHandler->getDPMBase()->speciesHandler.getMixedObject(getSpecies(),P->getSpecies());
675  const Mdouble sumOfInteractionRadii = getSumOfInteractionRadii(P);
676  if (distanceSquared >= sumOfInteractionRadii * sumOfInteractionRadii) {
677  return nullptr;
678  }
679  BaseInteraction* const C = interactionHandler->getInteraction(P, this, timeStamp);
680  const Mdouble distance = std::sqrt(distanceSquared);
681  C->setNormal(branchVector / distance);
682  C->setOverlap(P->getRadius() + getRadius() - distance);
683  C->setDistance(distance);
684  C->setContactPoint(P->getPosition() - (P->getRadius() - 0.5 * C->getOverlap()) * C->getNormal());
685  // The above contact point was chosen historically. It is not clear whether
686  // the current implementation or the two lines below are correct or incorrect.
687  // If solid arguments for the lines below are presented, it can be rediscussed.
688  // If so, also BaseWall::getInteractionWith should be adapted.
689  //Mdouble ratio=P->getRadius()/(getRadius()+P->getRadius());
690  //C->setContactPoint(P->getPosition() - (P->getRadius() - ratio * C->getOverlap()) * C->getNormal());
691  return C;
692 }
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
Matrix< Scalar, Dynamic, Dynamic > C
Definition: bench_gemm.cpp:49
Stores information about interactions between two interactable objects; often particles but could be ...
Definition: BaseInteraction.h:39
Mdouble getSumOfInteractionRadii(const BaseParticle *particle) const
returns the sum of the radii plus the interactionDistance
Definition: BaseParticle.h:362
BaseInteraction * getInteraction(BaseInteractable *P, BaseInteractable *I, unsigned timeStamp)
Returns the Interaction between the BaseInteractable's P and I.
Definition: InteractionHandler.cc:126
Mdouble getLengthSquared() const
Calculates the squared length of this Vec3D: .
Definition: Vector.cc:164
Definition: matrices.h:74
double P
Uniform pressure.
Definition: TwenteMeshGluing.cpp:77

References InteractionHandler::getInteraction(), Vec3D::getLengthSquared(), BaseInteractable::getPosition(), getRadius(), getSumOfInteractionRadii(), Global_Physical_Variables::P, and sqrt().

Referenced by SphericalSuperQuadricCollision::actionsAfterTimeStep(), ContactDetectionRotatedSpheresTest::actionsAfterTimeStep(), DPMBase::computeInternalForce(), PeriodicBoundaryHandler::processLocalInteractionData(), PeriodicBoundaryHandler::processReceivedInteractionData(), and Domain::processReceivedInteractionData().

◆ getInvInertia()

MatrixSymmetric3D BaseParticle::getInvInertia ( ) const
inline

Returns the inverse of the particle's inertia tensor.

Returns
the inverse of the particle's inertia tensor
270  { return invInertia_; }

References invInertia_.

Referenced by integrateAfterForceComputation(), and integrateBeforeForceComputation().

◆ getInvMass()

Mdouble BaseParticle::getInvMass ( ) const
inlineoverridevirtual

Returns the inverse of the particle's mass.

Returns
the inverse of the particle's mass

Reimplemented from BaseInteractable.

277  { return invMass_; }

References invMass_.

Referenced by integrateAfterForceComputation(), ClumpParticle::integrateAfterForceComputation(), integrateBeforeForceComputation(), and ClumpParticle::integrateBeforeForceComputation().

◆ getKineticEnergy()

Mdouble BaseParticle::getKineticEnergy ( ) const
virtual

Calculates the particle's translational kinetic energy.

Returns
the particle's translational kinetic energy

Calculates the particle's kinetic energy

Returns
the particle's kinetic energy

Reimplemented in NonSphericalParticle, and ClumpParticle.

447 {
448  if (isFixed())
449  return 0.0;
450  else
451  return 0.5 * getMass() * getVelocity().getLengthSquared();
452  //
453 }
static Mdouble getLengthSquared(const Vec3D &a)
Calculates the squared length of a Vec3D: .
Definition: Kernel/Math/Vector.h:324

References Vec3D::getLengthSquared(), getMass(), BaseInteractable::getVelocity(), and isFixed().

Referenced by NonSphericalParticle::getKineticEnergy().

◆ getMass()

◆ getMaxInteractionRadius()

Mdouble BaseParticle::getMaxInteractionRadius ( ) const
inline

Returns the particle's interaction radius, which might be different from radius_ (e.g., when dealing with wet particles)

Calculates the interaction radius of the particle (when it comes to interaction with other particles), including the effect of a possible additional 'interaction distance' besides the 'normal' radius. The interaction radius differs from the radius_, for example, when dealing with wet particles (i.e. particles with an additional liquid layer, which is dealt with in the particle's species).

Returns
the particle's interaction radius for particle-particle interaction
345  {
346  return getRadius() + getSpecies()->getMaxInteractionDistance() * 0.5;
347  }
Mdouble getMaxInteractionDistance() const
returns the largest separation distance at which adhesive short-range forces can occur.
Definition: ParticleSpecies.h:92

References ParticleSpecies::getMaxInteractionDistance(), getRadius(), and BaseInteractable::getSpecies().

Referenced by InsertionBoundary::checkBoundaryBeforeTimeStep(), Mercury3Dclump::checkClumpForInteractionPeriodic(), ParticleHandler::checkExtrema(), DPMBase::checkParticleForInteractionLocalPeriodic(), Mercury2D::computeInternalForces(), Mercury3D::computeInternalForces(), LeesEdwardsBoundary::createHorizontalPeriodicParticle(), ShearBoxBoundary::createHorizontalPeriodicParticle(), AngledPeriodicBoundary::createPeriodicParticle(), CircularPeriodicBoundary::createPeriodicParticle(), ConstantMassFlowMaserBoundary::createPeriodicParticle(), PeriodicBoundary::createPeriodicParticle(), SubcriticalMaserBoundary::createPeriodicParticle(), SubcriticalMaserBoundaryTEST::createPeriodicParticle(), TimeDependentPeriodicBoundary::createPeriodicParticle(), LeesEdwardsBoundary::createVerticalPeriodicParticle(), ShearBoxBoundary::createVerticalPeriodicParticle(), SuperQuadricParticle::getInteractionWith(), ParticleHandler::getLargestInteractionRadiusLocal(), ParticleHandler::getSmallestInteractionRadiusLocal(), Mercury3D::hGridFindContactsWithTargetCell(), Mercury2D::hGridFindParticleContacts(), Mercury3D::hGridFindParticleContacts(), Mercury2D::hGridGetInteractingParticleList(), Mercury3D::hGridGetInteractingParticleList(), Mercury2D::hGridHasParticleContacts(), Mercury3D::hGridHasParticleContacts(), MercuryBase::hGridNeedsRebuilding(), HGrid::insertParticleToHgrid(), and RotatingDrum::setupInitialConditions().

◆ getMomentum()

Vec3D BaseParticle::getMomentum ( ) const
inline
312  { return getVelocity() / invMass_; }

References BaseInteractable::getVelocity(), and invMass_.

◆ getName()

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

Returns the name of the object.

Returns the name of the object; in this case 'BaseParticle'.

Returns
The object name.

Implements BaseObject.

Reimplemented in NonSphericalParticle, SuperQuadricParticle, SphericalParticle, and ClumpParticle.

333 {
334  return "BaseParticle";
335 }

◆ getNameVTK()

std::string BaseParticle::getNameVTK ( unsigned  i) const
virtual

Reimplemented in ClumpParticle.

827 {
828  return "";
829 }

Referenced by ParticleVtkWriter::writeExtraFields().

◆ getNumberOfFieldsVTK()

unsigned BaseParticle::getNumberOfFieldsVTK ( ) const
virtual

Reimplemented in ClumpParticle.

817 {
818  return 0;
819 }

◆ getParticleDimensions()

unsigned int BaseParticle::getParticleDimensions ( ) const

Returns the particle's dimensions (either 2 or 3).

Returns the amount of dimensions of the particle (2 or 3, basically)

Returns
the number of dimension of the particle
765 {
767 }
DPMBase * getDPMBase()
Gets the problem that is solved using this handler.
Definition: BaseHandler.h:733
unsigned int getParticleDimensions() const
Returns the particle dimensionality.
Definition: DPMBase.cc:1458

References BaseHandler< T >::getDPMBase(), getHandler(), and DPMBase::getParticleDimensions().

Referenced by computeMass(), SuperQuadricParticle::computeMass(), and getVolume().

◆ getPeriodicComplexity() [1/2]

const std::vector< int > & BaseParticle::getPeriodicComplexity ( )

Obtains the periodic communication complexity of the particle.

Todo:
resolve this hack
215 {
217  //hack: generally you'd add particles after declaring the boundaries
218  //but no official programming guildelines rules have been setup for that
219  //So incase that doesnt happen we need to resize this periodicComplexity
220  if (periodicComplexity_.empty())
221  {
222  const unsigned numberOfPeriodicBoundaries = getHandler()->getDPMBase()->periodicBoundaryHandler.getSize();
223  if (numberOfPeriodicBoundaries > 0)
224  {
225  periodicComplexity_.resize(numberOfPeriodicBoundaries, 0);
226  }
227  }
228  return periodicComplexity_;
229 }
unsigned int getSize() const
Gets the size of the particleHandler (including mpi and periodic particles)
Definition: BaseHandler.h:663
PeriodicBoundaryHandler periodicBoundaryHandler
Internal handler that deals with periodic boundaries, especially in a parallel build.
Definition: DPMBase.h:1463

References BaseHandler< T >::getDPMBase(), getHandler(), BaseHandler< T >::getSize(), DPMBase::periodicBoundaryHandler, and periodicComplexity_.

Referenced by PeriodicBoundaryHandler::findNewParticle(), SubcriticalMaserBoundaryTEST::modifyGhostAfterCreation(), PeriodicBoundaryHandler::processLocalGhostParticles(), PeriodicBoundaryHandler::shiftParticle(), PeriodicBoundaryHandler::updateMaserParticle(), PeriodicBoundaryHandler::updateParticles(), and PeriodicBoundaryHandler::updateParticleStatus().

◆ getPeriodicComplexity() [2/2]

int BaseParticle::getPeriodicComplexity ( int  index)

Gets the periodic communication complexity of a certain boundary.

Todo:
TW Marnix, this is indeed a hack; you should call a setter every time you add a value to the periodic boundary handler (this function takes 0.5% cpu time in the speedtest)
232 {
233  //hack: generally you'd add particles after declaring the boundaries
234  //but no official programming guildelines rules have been setup for that
235  //So incase that doesnt happen we need to resize this periodicComplexity
237  if (periodicComplexity_.empty())
238  {
239  const unsigned numberOfPeriodicBoundaries = getHandler()->getDPMBase()->periodicBoundaryHandler.getSize();
240  if (numberOfPeriodicBoundaries > 0)
241  {
242  periodicComplexity_.resize(numberOfPeriodicBoundaries, 0);
243  }
244  }
245 
246  return periodicComplexity_[index];
247 }

References BaseHandler< T >::getDPMBase(), getHandler(), BaseHandler< T >::getSize(), DPMBase::periodicBoundaryHandler, and periodicComplexity_.

◆ getPeriodicFromParticle()

◆ getPreviousPeriodicComplexity()

const std::vector< int > & BaseParticle::getPreviousPeriodicComplexity ( ) const

Sets the previous periodic communication complexity of the particle.

255 {
257 }

References previousPeriodicComplexity_.

Referenced by PeriodicBoundaryHandler::updateParticles(), and PeriodicBoundaryHandler::updateParticleStatus().

◆ getPreviousPosition()

const Vec3D& BaseParticle::getPreviousPosition ( ) const
inline

Returns the particle's position in the previous time step.

Returns the particle's position in the previous time step.

Returns
(reference to) the previous position of the particle
379  { return previousPosition_; }

References previousPosition_.

Referenced by getDisplacement2(), and PeriodicBoundaryHandler::updateParticles().

◆ getRadius()

Mdouble BaseParticle::getRadius ( ) const
inline

Returns the particle's radius.

Returns
the particle's radius
332  { return radius_; }

References radius_.

Referenced by RotatingDrumBidisperseInitialise::actionsAfterTimeStep(), ConstantMassFlowMaserBoundary::activateMaser(), SilbertPeriodic::add_flow_particles(), statistics_while_running< T >::auto_set_domain(), statistics_while_running< T >::auto_set_z(), FixedClusterInsertionBoundary::checkBoundaryBeforeTimeStep(), InsertionBoundary::checkBoundaryBeforeTimeStep(), RandomClusterInsertionBoundary::checkBoundaryBeforeTimeStep(), Funnel::cleanChute(), CGHandler::computeContactPoints(), ChuteWithPeriodicInflow::computeInternalForces(), computeMass(), Funnel::create_funnel(), LawinenBox::create_inflow_particle(), ChuteWithContraction::create_inflow_particle(), AngleOfRepose::create_inflow_particle(), FlowRule::create_inflow_particle(), SilbertPeriodic::create_inflow_particle(), SegregationWithHopper::create_inflow_particle(), Slide::create_rough_wall(), Chute::createBottom(), CurvyChute::createBottom(), SphericalIndenter::getIndenterHeight(), SuperQuadricParticle::getInteractionRadius(), getInteractionWith(), getMaxInteractionRadius(), MarbleRun::getParticleMass(), ParticleHandler::getPSD(), ParticleParticleCollision::getRelativeVelocity(), WallParticleCollision::getRelativeVelocity(), getSumOfInteractionRadii(), getSurfaceArea(), getWallInteractionRadius(), LawinenBox::LawinenBox(), SphericalIndenter::outputXBallsData(), ChuteWithPeriodicInflow::outputXBallsDataParticlee(), SinterPair::printTime(), FileReader::read(), Slide::set_Walls(), SphericalIndenter::setIndenterHeight(), CFDDEMCoupleTest::setupInitialConditions(), ForceLawsMPI2Test::setupInitialConditions(), MarbleRun::setupInitialConditions(), statistics_while_running< T >::setupInitialConditions(), StressStrainControl::setupInitialConditions(), ParticleCreation::setupInitialConditions(), ParticleParticleCollision::setupInitialConditions(), WallParticleCollision::setupInitialConditions(), UnionOfWalls::setupInitialConditions(), MovingWallTangential::setupInitialConditions(), ChuteBottom::setupInitialConditions(), Slide::Slide(), statistics_while_running< T >::statistics_while_running(), SingleParticle< SpeciesType >::writeEneTimeStep(), and DPMBase::writeFstatHeader().

◆ getRotationalEnergy()

Mdouble BaseParticle::getRotationalEnergy ( ) const
virtual

Calculates the particle's rotational kinetic energy.

Returns
the particle's rotational kinetic energy

Reimplemented in NonSphericalParticle, and ClumpParticle.

456 {
457  if (isFixed())
458  return 0.0;
459  else
461 }
MatrixSymmetric3D getInertia() const
Definition: BaseParticle.h:314
static Mdouble dot(const Vec3D &a, const Vec3D &b)
Calculates the dot product of two Vec3D: .
Definition: Vector.cc:56

References Vec3D::dot(), BaseInteractable::getAngularVelocity(), getInertia(), and isFixed().

Referenced by NonSphericalParticle::getRotationalEnergy().

◆ getSumOfInteractionRadii()

Mdouble BaseParticle::getSumOfInteractionRadii ( const BaseParticle particle) const
inline

returns the sum of the radii plus the interactionDistance

362  {
363  return getRadius() + particle->getRadius() + getInteractionDistance((const BaseInteractable*)particle);
364  }
Defines the basic properties that a interactable object can have.
Definition: BaseInteractable.h:34
Mdouble getInteractionDistance(const BaseInteractable *i) const
Returns the interactionDistance_ of the mixed species of this particle and the particle or wall i.
Definition: BaseParticle.h:352

References getInteractionDistance(), and getRadius().

Referenced by Domain::findNewMPIInteractions(), SuperQuadricParticle::getInitialGuessForContact(), getInteractionWith(), and isInContactWith().

◆ getSurfaceArea()

virtual Mdouble BaseParticle::getSurfaceArea ( ) const
inlinevirtual

Reimplemented in NonSphericalParticle.

309  { return 4.0*constants::pi*getRadius()*getRadius(); }

References getRadius(), and constants::pi.

◆ getTimeStamp()

unsigned BaseParticle::getTimeStamp ( ) const
356 {
357  return timeStamp_;
358 }

References timeStamp_.

◆ getTypeVTK()

std::string BaseParticle::getTypeVTK ( unsigned  i) const
virtual

Reimplemented in ClumpParticle.

822 {
823  return "";
824 }

Referenced by ParticleVtkWriter::writeExtraFields().

◆ getVolume()

Mdouble BaseParticle::getVolume ( ) const
virtual

Get Particle volume function, which required a reference to the Species vector. It returns the volume of the Particle.

Returns the volume of the BaseParticle, which is calculated using its number of dimensions and radius.

Returns
The actual volume of this BaseParticle.

Reimplemented in SuperQuadricParticle.

127 {
128  if (handler_ == nullptr)
129  {
130  logger(ERROR, "[BaseParticle::getVolume] no particle handler specified");
131  return 0;
132  }
133  switch (getParticleDimensions())
134  {
135  case 3:
136  return (4.0 / 3.0 * constants::pi * radius_ * radius_ * radius_);
137  case 2:
138  return (constants::pi * radius_ * radius_);
139  case 1:
140  return (2.0 * radius_);
141  default:
142  logger(ERROR, "[BaseParticle::getVolume] dimension of the particle is not set");
143  return 0;
144  }
145 }
@ ERROR

References ERROR, getParticleDimensions(), handler_, logger, constants::pi, and radius_.

Referenced by InsertionBoundary::checkBoundaryBeforeTimeStep(), ScaleCoupling< M, O >::computeProjectionMatrix(), and ExtremeOverlapVolumeUnitTest::setupInitialConditions().

◆ getWallInteractionRadius()

Mdouble BaseParticle::getWallInteractionRadius ( const BaseWall wall) const
inline

returns the radius plus the interactionDistance

369  {
370  return getRadius() + getInteractionDistance((const BaseInteractable*)wall);
371  }

References getInteractionDistance(), and getRadius().

Referenced by TriangleMeshWall::getDistanceAndNormal().

◆ integrateAfterForceComputation()

void BaseParticle::integrateAfterForceComputation ( double  time,
double  timeStep 
)
virtual

Second step of Velocity Verlet integration.

Second step of Velocity Verlet integration (see also http://en.wikipedia.org/wiki/Verlet_integration#Velocity_Verlet).

Parameters
[in]timecurrent time
[in]timeStepcurrent time step

Reimplemented in ClumpParticle.

743 {
744  if (getInvMass() == 0.0)
745  {
746  //Updates a baseParticle with a prescribed motion
748  }
749  else
750  {
751  accelerate(getForce() * getInvMass() * 0.5 * timeStep);
752  if (getHandler()->getDPMBase()->getRotation())
753  {
755  getOrientation().rotateInverseInertiaTensor(getInvInertia()) * getTorque() * 0.5 * timeStep);
756  }
757  }
758 }
const Quaternion & getOrientation() const
Returns the orientation of this BaseInteractable.
Definition: BaseInteractable.h:209
const Vec3D & getForce() const
Returns the force on this BaseInteractable.
Definition: BaseInteractable.h:105
void integrateAfterForceComputation(double time, double timeStep)
This is part of the integration routine for objects with infinite mass.
Definition: BaseInteractable.cc:589
const Vec3D & getTorque() const
Returns the torque on this BaseInteractable.
Definition: BaseInteractable.h:117
MatrixSymmetric3D getInvInertia() const
Returns the inverse of the particle's inertia tensor.
Definition: BaseParticle.h:269
void angularAccelerate(const Vec3D &angVel)
Increases the particle's angularVelocity_ by the given vector.
Definition: BaseParticle.cc:630
void accelerate(const Vec3D &vel)
Increases the particle's velocity_ by the given vector.
Definition: BaseParticle.cc:620
Mdouble getInvMass() const override
Returns the inverse of the particle's mass.
Definition: BaseParticle.h:276

References accelerate(), angularAccelerate(), BaseInteractable::getForce(), getHandler(), getInvInertia(), getInvMass(), BaseInteractable::getOrientation(), BaseInteractable::getTorque(), and BaseInteractable::integrateAfterForceComputation().

◆ integrateBeforeForceComputation()

void BaseParticle::integrateBeforeForceComputation ( double  time,
double  timeStep 
)
virtual

First step of Velocity Verlet integration.

First step of Velocity Verlet integration (see also http://en.wikipedia.org/wiki/Verlet_integration#Velocity_Verlet).

Parameters
[in]timecurrent time
[in]timeStepcurrent time step
Todo:
If the position is described by the user, one should also call BaseInteractable::integrateBeforeForceComputation. To check if it works correctly, remove the p0.fixParticle() line from the DrivenParticleUnitTest
Author
irana

Reimplemented in ClumpParticle.

701 {
706 
707  if (getInvMass() == 0.0)
708  {
710  }
711  else
712  {
713  // For periodic particles in parallel the previous position is required
714  // also required by FluxBoundary
716 
717  accelerate(getForce() * getInvMass() * 0.5 * timeStep);
718  const Vec3D displacement = getVelocity() * timeStep;
719  move(displacement);
720  DPMBase* const dpm = getHandler()->getDPMBase();
721  if (!dpm->getHGridUpdateEachTimeStep())
722  {
723  dpm->hGridUpdateMove(this, displacement.getLengthSquared());
724  }
725 
726  if (dpm->getRotation())
727  {
729  getOrientation().rotateInverseInertiaTensor(getInvInertia()) * getTorque() * 0.5 * timeStep);
730  //apply to rotation quaternion q: q = normalise(q + \tilde{C}\omega*timeStep) (see Wouter's notes)
731  rotate(getAngularVelocity() * timeStep);
732  }
733  }
734 }
void integrateBeforeForceComputation(double time, double timeStep)
This is part of integrate routine for objects with infinite mass.
Definition: BaseInteractable.cc:516
virtual void move(const Vec3D &move)
Moves this BaseInteractable by adding an amount to the position.
Definition: BaseInteractable.cc:193
virtual void rotate(const Vec3D &angularVelocityDt)
Rotates this BaseInteractable.
Definition: BaseInteractable.cc:208
void setPreviousPosition(const Vec3D &pos)
Sets the particle's position in the previous time step.
Definition: BaseParticle.cc:600
The DPMBase header includes quite a few header files, defining all the handlers, which are essential....
Definition: DPMBase.h:56
virtual bool getHGridUpdateEachTimeStep() const
Definition: DPMBase.cc:1707
virtual void hGridUpdateMove(BaseParticle *, Mdouble)
Definition: DPMBase.cc:1922
bool getRotation() const
Indicates whether particle rotation is enabled or disabled.
Definition: DPMBase.h:554

References accelerate(), angularAccelerate(), BaseInteractable::getAngularVelocity(), BaseHandler< T >::getDPMBase(), BaseInteractable::getForce(), getHandler(), DPMBase::getHGridUpdateEachTimeStep(), getInvInertia(), getInvMass(), Vec3D::getLengthSquared(), BaseInteractable::getOrientation(), BaseInteractable::getPosition(), DPMBase::getRotation(), BaseInteractable::getTorque(), BaseInteractable::getVelocity(), DPMBase::hGridUpdateMove(), BaseInteractable::integrateBeforeForceComputation(), BaseInteractable::move(), BaseInteractable::rotate(), and setPreviousPosition().

◆ isClump()

◆ isFixed()

◆ isInContactWith()

bool BaseParticle::isInContactWith ( const BaseParticle P) const
virtual

Get whether or not this particle is in contact with the given particle.

Reimplemented in SuperQuadricParticle.

846 {
847  if (P->getName() != "Superquadric")
848  {
849  return Vec3D::getDistanceSquared(getPosition(), P->getPosition()) <
851  }
852  return P->isInContactWith(this);
853 }
static Mdouble getDistanceSquared(const Vec3D &a, const Vec3D &b)
Calculates the squared distance between two Vec3D: .
Definition: Kernel/Math/Vector.h:303

References Vec3D::getDistanceSquared(), BaseInteractable::getPosition(), getSumOfInteractionRadii(), Global_Physical_Variables::P, and mathsFunc::square().

Referenced by DPMBase::areInContact().

◆ isInMPIDomain()

bool BaseParticle::isInMPIDomain ( )

Indicates if the particle is in the communication zone of the mpi domain.

260 {
261  return isInMPIDomain_;
262 }

References isInMPIDomain_.

Referenced by Domain::findNewMPIInteractions(), Domain::findNewMPIParticle(), and PeriodicBoundaryHandler::updateParticleStatus().

◆ isInPeriodicDomain()

bool BaseParticle::isInPeriodicDomain ( ) const

Indicates if the particle is in the periodic boundary communication zone.

271 {
272  return isInPeriodicDomain_;
273 }

References isInPeriodicDomain_.

Referenced by PeriodicBoundaryHandler::checkIfAddNewParticle().

◆ isMaserParticle()

bool BaseParticle::isMaserParticle ( ) const

◆ isMPIParticle()

bool BaseParticle::isMPIParticle ( ) const

Indicates if this particle is a ghost in the MPI domain.

165 {
166  //make mpi-dependent so the compiler can optimise
167 #ifdef MERCURYDPM_USE_MPI
168  return isMPIParticle_;
169 #else
170  return false;
171 #endif
172 }

References isMPIParticle_.

Referenced by InteractionHandler::getLiquidBridgeVolume(), InteractionHandler::getNumberOfLiquidBridges(), DPMBase::outputXBallsData(), ParticleVtkWriter::particleMustBeWritten(), and PeriodicBoundaryHandler::updateParticleStatus().

◆ isPebble()

◆ isPeriodicGhostParticle()

bool BaseParticle::isPeriodicGhostParticle ( ) const

Indicates if this particle is a ghost in the periodic boundary.

281 {
283 }

References isPeriodicGhostParticle_.

Referenced by PeriodicBoundaryHandler::checkIfAddNewParticle(), Domain::findNewMPIParticle(), DPMBase::outputXBallsData(), and ParticleVtkWriter::particleMustBeWritten().

◆ isSphericalParticle()

virtual bool BaseParticle::isSphericalParticle ( ) const
pure virtual

this flag is used to decide whether to compute orientation, which is not necessary for spherical particles

Todo:
This flag is used badly, and used to determine whether particles are superquadric

Implemented in SphericalParticle, NonSphericalParticle, and ClumpParticle.

◆ movePrevious()

void BaseParticle::movePrevious ( const Vec3D posMove)

Adds a vector to the particle's previousPosition_.

Lets you add a vector to the particle's previousPosition_ vector.

Parameters
[in]posMovethe vector to be added to the current previousPosition_ vector.
611 {
612  previousPosition_ += posMove;
613 }

References previousPosition_.

◆ oldRead()

void BaseParticle::oldRead ( std::istream &  is)
virtual

Should NOT BE USED by any user, only used to read old restart files! Is expected to be obsolete by Mercury 2.0. Please use BaseParticle::read() instead.

This is the previously used version of the read function. Now just kept for legacy purposes.

Deprecated:
Should be gone in Mercury 2.0. Use BaseParticle::read() instead.
Todo:
incorporate contact information
382 {
383  logger(DEBUG, "reading particle old-style");
384  static unsigned int id = 0;
385  unsigned int indSpecies = 0;
386  unsigned int numberOfContacts = 0;
387  Vec3D orientation;
388  Vec3D position;
389  Vec3D velocity;
390  Vec3D angularVelocity;
391  double invInertiaScalar;
392  double dummy = 0;
393  is >> position >> velocity >> radius_ >> orientation >> angularVelocity;
394  is >> invMass_ >> invInertiaScalar >> numberOfContacts;
396  for (unsigned int i = 0; i < 12 * numberOfContacts; ++i)
397  {
398  is >> dummy;
399  }
400  is >> indSpecies;
401  setPosition(position);
403  Quaternion q;
404  q.setEuler(orientation);
405  setOrientation(q);
406  setAngularVelocity(angularVelocity);
407  invInertia_.XX = invInertiaScalar;
408  invInertia_.YY = invInertiaScalar;
409  invInertia_.ZZ = invInertiaScalar;
410 
412  setId(id);
413  setIndex(id);
414  id++;
415 }
virtual void setOrientation(const Quaternion &orientation)
Sets the orientation of this BaseInteractable.
Definition: BaseInteractable.h:239
virtual void setPosition(const Vec3D &position)
Sets the position of this BaseInteractable.
Definition: BaseInteractable.h:218
virtual void setIndSpecies(unsigned int indSpecies)
Sets the index of the Species of this BaseInteractable.
Definition: BaseInteractable.h:77
void setIndex(unsigned int index)
Allows one to assign an index to an object in the handler/container.
Definition: BaseObject.cc:42
void setId(unsigned long id)
Assigns a unique identifier to each object in the handler (container) which remains constant even aft...
Definition: BaseObject.cc:50
Mdouble ZZ
Definition: MatrixSymmetric.h:21
Mdouble YY
Definition: MatrixSymmetric.h:21
Mdouble XX
The six distinctive matrix elements.
Definition: MatrixSymmetric.h:21
This class contains the 4 components of a quaternion and the standard operators and functions needed ...
Definition: Kernel/Math/Quaternion.h:42
EIGEN_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:2019
double velocity(const double &t)
Angular velocity as function of time t.
Definition: jeffery_orbit.cc:107

References DEBUG, i, invInertia_, invMass_, logger, Eigen::numext::q, radius_, BaseInteractable::setAngularVelocity(), BaseObject::setId(), BaseObject::setIndex(), BaseInteractable::setIndSpecies(), BaseInteractable::setOrientation(), BaseInteractable::setPosition(), BaseInteractable::setVelocity(), Jeffery_Solution::velocity(), MatrixSymmetric3D::XX, MatrixSymmetric3D::YY, and MatrixSymmetric3D::ZZ.

◆ printHGrid()

void BaseParticle::printHGrid ( std::ostream &  os) const

Adds particle's HGrid level and cell coordinates to an ostream.

Adds the particle's HGridLevel_ and HGRid x/y/z positions to an std::ostream.

Parameters
[in,out]osthe ostream which has the mentioned properties added.
423 {
424  os << "Particle( HGRID_Level:" << hGridCell.getHGridLevel()
425  << ", HGRID_x:" << hGridCell.getHGridX()
426  << ", HGRID_y:" << hGridCell.getHGridY()
427  << ", HGRID_z:" << hGridCell.getHGridZ()
428  << ")";
429 }

References HGridCell::getHGridLevel(), HGridCell::getHGridX(), HGridCell::getHGridY(), HGridCell::getHGridZ(), and hGridCell.

◆ read()

void BaseParticle::read ( std::istream &  is)
overridevirtual

Particle read function, which accepts an std::istream as input.

Particle read function. Has an std::istream as argument, from which it extracts the radius_, invMass_ and invInertia_, respectively. From these the mass and inertia are deduced. An additional set of properties is read through the call to the parent's method BaseInteractable::read().

Parameters
[in,out]isinput stream with particle properties.

Reimplemented from BaseInteractable.

Reimplemented in SuperQuadricParticle, and ClumpParticle.

369 {
371  std::string dummy;
372  is >> dummy >> radius_ >> dummy >> invMass_;// >> dummy >> invInertia_;
373  helpers::readOptionalVariable(is, "timeStamp", timeStamp_); // Optional, for backwards compatibility.
374 }
void read(std::istream &is) override
Reads a BaseInteractable from an input stream.
Definition: BaseInteractable.cc:222
bool readOptionalVariable(std::istream &is, const std::string &name, T &variable)
Reads optional variables in the restart file.
Definition: FileIOHelpers.h:68
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286

References invMass_, radius_, BaseInteractable::read(), helpers::readOptionalVariable(), oomph::Global_string_for_annotation::string(), and timeStamp_.

Referenced by ClumpParticle::read(), SuperQuadricParticle::read(), and ParticleHandler::readAndCreateObject().

◆ setAxes()

virtual void BaseParticle::setAxes ( const Vec3D axes)
inlinevirtual

Only ustilised in case of superquadric particles.

Reimplemented in SuperQuadricParticle.

502 { }

◆ setCommunicationComplexity()

void BaseParticle::setCommunicationComplexity ( unsigned  complexity)

Set the communication complexity of the particle.

180 {
181  communicationComplexity_ = complexity;
182 }

References communicationComplexity_.

Referenced by Domain::addParticlesToLists(), and Domain::updateParticles().

◆ setExponents()

virtual void BaseParticle::setExponents ( const Mdouble eps1,
const Mdouble eps2 
)
inlinevirtual

Only ustilised in case of superquadric particles.

Reimplemented in SuperQuadricParticle.

507 {}

◆ setHandler()

void BaseParticle::setHandler ( ParticleHandler handler)

Sets the pointer to the particle's ParticleHandler.

Assigns the particle to a ParticleHandler, and assigns a species to it based on the particles indSpecies_ (BaseInteractable data member).

Parameters
[in]handlerpointer to the ParticleHandler
641 {
642  handler_ = handler;
643  setSpecies(getHandler()->getDPMBase()->speciesHandler.getObject(getIndSpecies()));
644 }
T * getObject(const unsigned int id)
Gets a pointer to the Object at the specified index in the BaseHandler.
Definition: BaseHandler.h:621
unsigned int getIndSpecies() const
Returns the index of the species associated with the interactable object.
Definition: BaseInteractable.h:67

References getHandler(), BaseInteractable::getIndSpecies(), BaseHandler< T >::getObject(), handler_, and setSpecies().

Referenced by InsertionBoundary::checkBoundaryBeforeTimeStep(), RandomClusterInsertionBoundary::checkBoundaryBeforeTimeStep(), Chute::createBottom(), ParticleHandler::readAndCreateObject(), setSpecies(), SphericalIndenter::setupInitialConditions(), and ChuteBottom::setupInitialConditions().

◆ setHGridLevel()

void BaseParticle::setHGridLevel ( const unsigned int  level)
inline

Sets the particle's HGrid level.

Parameters
[in]levelthe particle's HGrid level
448  { hGridCell.setHGridLevel(level); }

References hGridCell, and HGridCell::setHGridLevel().

Referenced by HGrid::insertParticleToHgrid().

◆ setHGridNextObject()

void BaseParticle::setHGridNextObject ( BaseParticle p)
inline

Sets the pointer to the next object in the particle's HGrid cell & level.

Parameters
[in]ppointer to the next object
456  { hGridNextObject_ = p; }

References hGridNextObject_, and p.

Referenced by Mercury2D::hGridRemoveParticle(), Mercury3D::hGridRemoveParticle(), Mercury2D::hGridUpdateParticle(), and Mercury3D::hGridUpdateParticle().

◆ setHGridPrevObject()

void BaseParticle::setHGridPrevObject ( BaseParticle p)
inline

Sets the pointer to the previous object in the particle's HGrid cell & level.

Parameters
[in]ppointer to the previous object
464  { hGridPrevObject_ = p; }

References hGridPrevObject_, and p.

Referenced by Mercury2D::hGridRemoveParticle(), Mercury3D::hGridRemoveParticle(), Mercury2D::hGridUpdateParticle(), and Mercury3D::hGridUpdateParticle().

◆ setHGridX()

void BaseParticle::setHGridX ( const int  x)
inline

Sets the particle's HGrid cell X-coordinate.

Set the x-index of the particle's hGrid cell position

Parameters
[in]xx-index of particle's HGrid cell
425  { hGridCell.setHGridX(x); }
list x
Definition: plotDoE.py:28

References hGridCell, HGridCell::setHGridX(), and plotDoE::x.

Referenced by Mercury2D::hGridUpdateParticle(), and Mercury3D::hGridUpdateParticle().

◆ setHGridY()

void BaseParticle::setHGridY ( const int  y)
inline

Sets the particle's HGrid cell Y-coordinate.

Set the y-index of the particle's hGrid cell position

Parameters
[in]yy-index of particle's HGrid cell
433  { hGridCell.setHGridY(y); }
Scalar * y
Definition: level1_cplx_impl.h:128

References hGridCell, HGridCell::setHGridY(), and y.

Referenced by Mercury2D::hGridUpdateParticle(), and Mercury3D::hGridUpdateParticle().

◆ setHGridZ()

void BaseParticle::setHGridZ ( const int  z)
inline

Sets the particle's HGrid cell Z-coordinate.

Set the y-index of the particle's hGrid cell position

Parameters
[in]zz-index of particle's HGrid cell
441  { hGridCell.setHGridZ(z); }

References hGridCell, and HGridCell::setHGridZ().

Referenced by Mercury3D::hGridUpdateParticle().

◆ setIndSpecies()

void BaseParticle::setIndSpecies ( unsigned int  indSpecies)
overridevirtual
Deprecated:
Please use setSpecies(const ParticleSpecies*) instead.
Todo:
MX: this index is used in the MPI transmission. This should be "undeprecated"

Set the particle's species and species' index. Logs a warning if no ParticleHandler is assigned.

Parameters
[in]indSpeciesThe index of the species in the SpeciesHandler.
Todo:
TW do we have to update the species stored in the interactions here?

Reimplemented from BaseInteractable.

775 {
776  if (handler_ != nullptr)
777  {
778  //BaseInteractable::setIndSpecies(indSpecies);
781  }
782  else
783  {
785  logger(ERROR, "setIndSpecies called on a particle with no particle handler.\n"
786  "Therefore I can't request the given species from the species handler.\n"
787  " PartID = %", getId());
788  }
789 }
SpeciesHandler speciesHandler
A handler to that stores the species type i.e. LinearViscoelasticSpecies, etc.
Definition: DPMBase.h:1433

References ERROR, BaseHandler< T >::getDPMBase(), BaseObject::getId(), BaseHandler< T >::getObject(), handler_, logger, BaseInteractable::setIndSpecies(), setSpecies(), and DPMBase::speciesHandler.

◆ setInertia() [1/2]

void BaseParticle::setInertia ( )
virtual

◆ setInertia() [2/2]

void BaseParticle::setInertia ( MatrixSymmetric3D  inertia)

Sets the particle's inertia_ (and adjusts invInertia_ accordingly)

Sets the particle's inertia and invInertia_.

Parameters
[in]newInertiathe new inertia to be set.
510 {
512 }
MatrixSymmetric3D inverse() const
Computes the inverse of a matrix; exits if the inverse doesn't exist.
Definition: MatrixSymmetric.cc:277

References MatrixSymmetric3D::inverse(), and invInertia_.

◆ setInfiniteInertia()

void BaseParticle::setInfiniteInertia ( )

Sets the particle's inertia_ to 'infinite' (1e20) and its invInertia_ to 0.

Sets the inertia to 1e20 and the invInertia_ (which is actually used in the calculations) to 0.

524 {
526 
527 // Clump* f;
528 // f->invInertia2_.setZero();
529 } //> i.e. no rotations
void setZero()
Sets all elements to zero.
Definition: MatrixSymmetric.cc:49

References invInertia_, and MatrixSymmetric3D::setZero().

Referenced by Slide::actionsBeforeTimeLoop().

◆ setInfo()

void BaseParticle::setInfo ( Mdouble  info)
virtual

Sets some user-defined information about this object (by default, species ID).

338 {
339  info_ = info;
340 }
int info
Definition: level2_cplx_impl.h:39

References info, and info_.

◆ setInMPIDomain()

void BaseParticle::setInMPIDomain ( bool  flag)

Flags the status of the particle if wether it is in the communication zone or not.

265 {
266  isInMPIDomain_ = flag;
267 }

References isInMPIDomain_.

Referenced by PeriodicBoundaryHandler::setMPIFlags(), Domain::updateParticles(), and PeriodicBoundaryHandler::updateParticleStatus().

◆ setInPeriodicDomain()

void BaseParticle::setInPeriodicDomain ( bool  flag)

◆ setInverseInertia()

void BaseParticle::setInverseInertia ( MatrixSymmetric3D  inverseInertia)

Sets the particle's inertia_ (and adjusts invInertia_ accordingly)

515 {
516  invInertia_ = inverseInertia;
517 }

References invInertia_.

◆ setMaserParticle()

void BaseParticle::setMaserParticle ( bool  flag)

Flags the status of the particle if it belongs to the maser boundary or not.

296 {
297  isMaserParticle_ = flag;
298 }

References isMaserParticle_.

Referenced by SubcriticalMaserBoundaryTEST::checkBoundaryAfterParticleMoved(), and PeriodicBoundaryHandler::updateMaserParticle().

◆ setMass()

void BaseParticle::setMass ( Mdouble  mass)

Sets the particle's mass.

Deprecated:
Please do not set the mass yourself, but use ParticleSpecies->computeMass instead. That makes sure

Sets the mass of the particle

Parameters
[in]massthe new particle's mass
564 {
565  logger(WARN, "WARNING: Do not use particle->setMass, instead use "
566  "particleSpecies->computeMass, since this function can cause "
567  "inconsistencies between the mass, density and radius of this particle!");
568  logger.assert_always(mass > 0.0 && !isFixed(),
569  "Error in BaseParticle::setMass, the given mass to be set must be positive.");
570 
571  invMass_ = 1.0 / mass;
572 }
@ WARN

References invMass_, isFixed(), logger, and WARN.

◆ setMassForP3Statistics()

void BaseParticle::setMassForP3Statistics ( Mdouble  mass)

Sets the particle's mass This function should not be used, but is necessary to extend the CG toolbox to non-spherical particles.

Sets the mass of the particle

Parameters
[in]massthe new particle's mass
585 {
586  if (mass > 0.0 && !isFixed())
587  {
588  invMass_ = 1.0 / mass;
589  }
590  else
591  {
592  logger(ERROR, "Error in BaseParticle::setMass, the given mass to be set must be positive.");
593  }
594 }

References ERROR, invMass_, isFixed(), and logger.

◆ setMPIParticle()

void BaseParticle::setMPIParticle ( bool  flag)

Flags the mpi particle status.

175 {
176  isMPIParticle_ = flag;
177 }

References isMPIParticle_.

Referenced by PeriodicBoundaryHandler::setMPIFlags(), and Domain::updateParticles().

◆ setPeriodicComplexity() [1/2]

void BaseParticle::setPeriodicComplexity ( int  index,
int  value 
)

Set the periodic communication complexity of the particle.

196 {
197  //hack: generally you'd add particles after declaring the boundaries
198  //but no official programming guildelines rules have been setup for that
199  //So incase that doesnt happen we need to resize this periodicComplexity
200  if (periodicComplexity_.empty())
201  {
202  int numberOfPeriodicBoundaries = getHandler()->getDPMBase()->periodicBoundaryHandler.getSize();
203  if (numberOfPeriodicBoundaries > 0)
204  {
205  //First initialisation of the periodic complexity assumes the particle is completely
206  //within the real domain
207  periodicComplexity_ = std::vector<int>(numberOfPeriodicBoundaries, 2);
208  }
209  }
210 
211  periodicComplexity_[index] = value;
212 }
squared absolute value
Definition: GlobalFunctions.h:87

References BaseHandler< T >::getDPMBase(), getHandler(), BaseHandler< T >::getSize(), DPMBase::periodicBoundaryHandler, periodicComplexity_, and Eigen::value.

◆ setPeriodicComplexity() [2/2]

void BaseParticle::setPeriodicComplexity ( std::vector< int complexity)

◆ setPeriodicFromParticle()

◆ setPeriodicGhostParticle()

void BaseParticle::setPeriodicGhostParticle ( bool  flag)

◆ setPreviousPeriodicComplexity()

void BaseParticle::setPreviousPeriodicComplexity ( std::vector< int complexity)

Set the previous periodic communication complexity of the paritcle.

250 {
251  previousPeriodicComplexity_ = complexity;
252 }

References previousPeriodicComplexity_.

Referenced by PeriodicBoundaryHandler::updateParticles().

◆ setPreviousPosition()

void BaseParticle::setPreviousPosition ( const Vec3D pos)

Sets the particle's position in the previous time step.

This is used to set the particle's previous position

Parameters
[in]posthe particle's previous position vector.
601 {
602  previousPosition_ = pos;
603 }

References previousPosition_.

Referenced by CircularPeriodicBoundary::checkBoundaryAfterParticleMoved(), integrateBeforeForceComputation(), ClumpParticle::integrateBeforeForceComputation(), Domain::processReceivedBoundaryParticleData(), and PeriodicBoundaryHandler::updateParticles().

◆ setRadius()

void BaseParticle::setRadius ( Mdouble  radius)
virtual

Sets the particle's radius_ (and adjusts the mass_ accordingly, based on the particle's species)

Sets the radius of the particle, and from that computes the new mass (using its species) and checks whether it now is either the smallest or biggest particle in its ParticleHandler.

Parameters
[in]radiusthe new radius

Reimplemented in SuperQuadricParticle.

549 {
550  radius_ = radius;
551  if (getHandler())
552  {
553  getSpecies()->computeMass(this);
554  getHandler()->checkExtrema(this);
555  }
556 
557 }
void checkExtrema(BaseParticle *P)
Checks if the extrema of this ParticleHandler needs updating.
Definition: ParticleHandler.cc:1171
void computeMass(BaseParticle *p) const
Compute Particle mass function, which required a reference to the Species vector. It computes the Par...
Definition: ParticleSpecies.cc:147
radius
Definition: UniformPSDSelfTest.py:15

References ParticleHandler::checkExtrema(), ParticleSpecies::computeMass(), getHandler(), BaseInteractable::getSpecies(), UniformPSDSelfTest::radius, and radius_.

Referenced by ClumpParticle::ClumpParticle(), Funnel::create_funnel(), LawinenBox::create_inflow_particle(), ChuteWithContraction::create_inflow_particle(), Funnel::create_inflow_particle(), AngleOfRepose::create_inflow_particle(), FlowRule::create_inflow_particle(), SilbertPeriodic::create_inflow_particle(), SegregationWithHopper::create_inflow_particle(), Chute::createBottom(), CurvyChute::createBottom(), BasicIntersectionOfWalls::getDistanceAndNormal(), BasicUnionOfWalls::getDistanceAndNormal(), TriangleMeshWall::getDistanceAndNormal(), BasicIntersectionOfWalls::getVTK(), BasicUnionOfWalls::getVTK(), main(), regimeForceUnitTest::regimeForceUnitTest(), RotatingDrumBidisperseInitialise::RotatingDrumBidisperseInitialise(), SuperQuadricParticle::setBoundingRadius(), MarbleRun::setParticleRadius(), RandomClusterInsertionBoundarySelfTest::setupInitialConditions(), CFDDEMCoupleTest::setupInitialConditions(), my_problem::setupInitialConditions(), ForceLawsMPI2Test::setupInitialConditions(), TwoByTwoMPIDomainMPI4Test::setupInitialConditions(), ParameterStudy1DDemo::setupInitialConditions(), ParameterStudy2DDemo::setupInitialConditions(), ParameterStudy3DDemo::setupInitialConditions(), StressStrainControl::setupInitialConditions(), Cstatic2d::setupInitialConditions(), Cstatic3D::setupInitialConditions(), GetDistanceAndNormalForIntersectionOfWalls::setupInitialConditions(), GetDistanceAndNormalForScrew::setupInitialConditions(), GetDistanceAndNormalForTriangleWall::setupInitialConditions(), Drum::setupInitialConditions(), HertzSelfTest::setupInitialConditions(), MindlinSelfTest::setupInitialConditions(), Penetration::setupInitialConditions(), Silo::setupInitialConditions(), ConstantMassFlowMaserSelfTest::setupInitialConditions(), WallParticleCollision::setupInitialConditions(), GetDistanceAndNormalForTriangleWalls::setupInitialConditions(), RollingOverTriangleWalls::setupInitialConditions(), UnionOfWalls::setupInitialConditions(), SphericalIndenter::setupInitialConditions(), SlidingSpheresUnitTest::setupInitialConditions(), ChargedBondedParticleUnitTest::setupInitialConditions(), HertzianSinterForceUnitTest::setupInitialConditions(), MovingIntersectionOfWallsUnitTest_Basic::setupInitialConditions(), MovingWall::setupInitialConditions(), MultiParticlesInsertion::setupInitialConditions(), PeriodicWallsWithSlidingFrictionUnitTest::setupInitialConditions(), PlasticForceUnitTest::setupInitialConditions(), my_problem_HGRID::setupInitialConditions(), SinterForceUnitTest::setupInitialConditions(), SpeciesTest::setupInitialConditions(), TangentialSpringUnitTest::setupInitialConditions(), ChuteBottom::setupInitialConditions(), and viscoElasticUnitTest::viscoElasticUnitTest().

◆ setSpecies()

void BaseParticle::setSpecies ( const ParticleSpecies species)
virtual

In addition to the functionality of BaseInteractable::setSpecies, this function sets the pointer to the particleHandler, which is needed to retrieve species information.

Todo:
TW: this function should also check if the particle is the correct particle for the species type

Sets the particle's species. If this particle does not have a handler yet, this function also assigns the ParticleHandler in the same DPMBase as the SpeciesHandler of the given species as its handler.

Parameters
[in]speciespointer to the ParticleSpecies object, to be set as the particle's species.
Todo:
TW should we chaeck here if we have the right kind of species for the right kind of particle?
Todo:
maybe these if statements should throw warnings
799 {
802  //set pointer to the ParticleHandler handler_, which is needed to retrieve
803  //species information
805  if (handler_ == nullptr)
806  {
807  SpeciesHandler* sH = species->getHandler();
808  DPMBase* dB = sH->getDPMBase();
809  if (dB != nullptr)
810  {
812  }
813  }
814 }
void setSpecies(const ParticleSpecies *species)
Sets the species of this BaseInteractable.
Definition: BaseInteractable.cc:163
void setHandler(ParticleHandler *handler)
Sets the pointer to the particle's ParticleHandler.
Definition: BaseParticle.cc:640
SpeciesHandler * getHandler() const
Returns the pointer to the handler to which this species belongs.
Definition: BaseSpecies.cc:78
ParticleHandler particleHandler
An object of the class ParticleHandler, contains the pointers to all the particles created.
Definition: DPMBase.h:1443
Container to store all ParticleSpecies.
Definition: SpeciesHandler.h:15

References BaseHandler< T >::getDPMBase(), BaseSpecies::getHandler(), handler_, DPMBase::particleHandler, setHandler(), and BaseInteractable::setSpecies().

Referenced by ChuteWithPeriodicInflow::AddContinuingBottom(), BaseParticle(), BoundariesSelfTest::BoundariesSelfTest(), ConstantMassFlowMaserBoundary::checkBoundaryAfterParticleMoved(), SubcriticalMaserBoundary::checkBoundaryAfterParticleMoved(), ChuteWithPeriodicInflowAndContinuingBottom::ChuteWithPeriodicInflowAndContinuingBottom(), ChuteWithPeriodicInflowAndContraction::ChuteWithPeriodicInflowAndContraction(), ChuteWithPeriodicInflowAndVariableBottom::ChuteWithPeriodicInflowAndVariableBottom(), ContractionWithPeriodicInflow::ContractionWithPeriodicInflow(), FlowRule::create_inflow_particle(), SilbertPeriodic::create_inflow_particle(), AngleOfRepose::createBaseSpecies(), SilbertPeriodic::createBaseSpecies(), Chute::createBottom(), CurvyChute::createBottom(), FluxBoundarySelfTest::FluxBoundarySelfTest(), BasicIntersectionOfWalls::getDistanceAndNormal(), BasicUnionOfWalls::getDistanceAndNormal(), TriangleMeshWall::getDistanceAndNormal(), BasicIntersectionOfWalls::getVTK(), BasicUnionOfWalls::getVTK(), Hertzian2DUnitTest::Hertzian2DUnitTest(), ChuteWithPeriodicInflow::integrateBeforeForceComputation(), LawinenBox::LawinenBox(), main(), MercuryCGSelfTest::MercuryCGSelfTest(), protectiveWall::protectiveWall(), InsertionBoundary::read(), ParticleHandler::readAndCreateObject(), regimeForceUnitTest::regimeForceUnitTest(), RotatingDrumBidisperseInitialise::RotatingDrumBidisperseInitialise(), setHandler(), setIndSpecies(), MeltableParticle::setSpecies(), CFDDEMCoupleTest::setupInitialConditions(), ForceLawsMPI2Test::setupInitialConditions(), TwoByTwoMPIDomainMPI4Test::setupInitialConditions(), SilbertPeriodic::setupInitialConditions(), GetDistanceAndNormalForIntersectionOfWalls::setupInitialConditions(), GetDistanceAndNormalForScrew::setupInitialConditions(), GetDistanceAndNormalForTriangleWall::setupInitialConditions(), Drum::setupInitialConditions(), HertzSelfTest::setupInitialConditions(), MindlinSelfTest::setupInitialConditions(), Penetration::setupInitialConditions(), Silo::setupInitialConditions(), ConstantMassFlowMaserSelfTest::setupInitialConditions(), InsertionBoundarySelfTest::setupInitialConditions(), GetDistanceAndNormalForTriangleWalls::setupInitialConditions(), RollingOverTriangleWalls::setupInitialConditions(), UnionOfWalls::setupInitialConditions(), SphericalIndenter::setupInitialConditions(), GranularCollapse::setupInitialConditions(), SlidingSpheresUnitTest::setupInitialConditions(), ChargedBondedParticleUnitTest::setupInitialConditions(), RandomClusterInsertionBoundarySelfTest::setupInitialConditions(), FullRestartTest::setupInitialConditions(), HertzianSinterForceUnitTest::setupInitialConditions(), MultiParticlesInsertion::setupInitialConditions(), PeriodicWallsWithSlidingFrictionUnitTest::setupInitialConditions(), PlasticForceUnitTest::setupInitialConditions(), SinterForceUnitTest::setupInitialConditions(), SpeciesTest::setupInitialConditions(), TangentialSpringUnitTest::setupInitialConditions(), Chute::setupInitialConditions(), ChuteBottom::setupInitialConditions(), ContactDetectionWithWallTester::setupParticleAndWall(), ContactDetectionTester::setupParticles(), SilbertPeriodic::SilbertPeriodic(), T_protectiveWall::T_protectiveWall(), and viscoElasticUnitTest::viscoElasticUnitTest().

◆ setTimeStamp()

void BaseParticle::setTimeStamp ( unsigned  timeStamp)
351 {
352  timeStamp_ = timeStamp;
353 }

References timeStamp_.

◆ unfix()

void BaseParticle::unfix ( )

Unfix Particle function, which required a reference to the Species vector. It unfixes a Particle by computing the Particles mass and inertia.

Unfixes the particle by computing the Particles mass and inertia, using the species and radius.

306 {
307  invMass_ = 1.0;
308  getSpecies()->computeMass(this);
309  if (getHandler())
311 }

References ParticleSpecies::computeMass(), getHandler(), BaseInteractable::getSpecies(), invMass_, and ParticleHandler::removedFixedParticle().

◆ write()

void BaseParticle::write ( std::ostream &  os) const
overridevirtual

Particle print function, which accepts an std::ostream as input.

BaseParticle print method, which accepts an os std::ostream as input. It prints human readable BaseParticle information to the std::ostream.

Parameters
[in,out]osstream to which the info is written

Reimplemented from BaseInteractable.

Reimplemented in SuperQuadricParticle, and ClumpParticle.

320 {
322  os << " radius " << radius_
323  << " invMass " << invMass_
324  << " timeStamp " << timeStamp_;
325  //invMass_ is a computed value, but needs to be stored to see if a particle is fixed
326 }
void write(std::ostream &os) const override
Write a BaseInteractable to an output stream.
Definition: BaseInteractable.cc:252

References invMass_, radius_, timeStamp_, and BaseInteractable::write().

Referenced by ClumpParticle::write(), and SuperQuadricParticle::write().

Friends And Related Function Documentation

◆ ParticleSpecies::computeMass

void ParticleSpecies::computeMass ( BaseParticle ) const
friend

Particle's position at previous time step.

Since ParticleSpecies is allowed to set the mass of a BaseParticle, it is a friend of this class.

Member Data Documentation

◆ clumpParticle_

BaseParticle* BaseParticle::clumpParticle_

Function that updates necessary quantities of a clump particle after adding a pebble.

Referenced by ClumpParticle::ClumpParticle(), getClump(), and ClumpParticle::setClump().

◆ communicationComplexity_

unsigned BaseParticle::communicationComplexity_
private

returns true if it flagged as being in MPI domain

Referenced by BaseParticle(), getCommunicationComplexity(), and setCommunicationComplexity().

◆ handler_

ParticleHandler* BaseParticle::handler_
private

Inverse Particle inverse inertia (for computation optimization)

Pointer to the particle's ParticleHandler container

Referenced by BaseParticle(), getHandler(), getVolume(), setHandler(), setIndSpecies(), and setSpecies().

◆ hGridCell

HGridCell BaseParticle::hGridCell
private

All hGrid-information: the cell (x,y,z,level), and the previous and next particle in this cell compared to this particle

Referenced by BaseParticle(), getHGridCell(), getHGridLevel(), getHGridX(), getHGridY(), getHGridZ(), printHGrid(), setHGridLevel(), setHGridX(), setHGridY(), and setHGridZ().

◆ hGridNextObject_

BaseParticle* BaseParticle::hGridNextObject_
private

◆ hGridPrevObject_

BaseParticle* BaseParticle::hGridPrevObject_
private

Pointer to the next Particle in the same HGrid cell.

Referenced by BaseParticle(), getHGridPrevObject(), and setHGridPrevObject().

◆ info_

Mdouble BaseParticle::info_
private

Referenced by BaseParticle(), getInfo(), and setInfo().

◆ invInertia_

◆ invMass_

◆ isClump_

◆ isInMPIDomain_

bool BaseParticle::isInMPIDomain_
private

returns true if the particle acts as an MPI particle instead of a real particle

Referenced by BaseParticle(), isInMPIDomain(), and setInMPIDomain().

◆ isInPeriodicDomain_

bool BaseParticle::isInPeriodicDomain_
private

◆ isMaserParticle_

bool BaseParticle::isMaserParticle_
private

Indicates the periodic complexity at current time step. Used to update periodic status.

Referenced by BaseParticle(), isMaserParticle(), and setMaserParticle().

◆ isMPIParticle_

bool BaseParticle::isMPIParticle_
private

Pointer to originating Particle.

Referenced by BaseParticle(), isMPIParticle(), and setMPIParticle().

◆ isPebble_

bool BaseParticle::isPebble_

◆ isPeriodicGhostParticle_

bool BaseParticle::isPeriodicGhostParticle_
private

bool that indicates if a particle is in the periodic domain of any boundary

Referenced by BaseParticle(), isPeriodicGhostParticle(), and setPeriodicGhostParticle().

◆ periodicComplexity_

std::vector<int> BaseParticle::periodicComplexity_
private

Indicates the periodic complexity at previous time step.

Referenced by BaseParticle(), getPeriodicComplexity(), and setPeriodicComplexity().

◆ periodicFromParticle_

BaseParticle* BaseParticle::periodicFromParticle_
private

Pointer to the previous Particle in the same HGrid cell.

Particle attributes

Referenced by BaseParticle(), getPeriodicFromParticle(), and setPeriodicFromParticle().

◆ previousPeriodicComplexity_

std::vector<int> BaseParticle::previousPeriodicComplexity_
private

Indicates if the particle is a ghost particle of a periodic particle.

Referenced by BaseParticle(), getPreviousPeriodicComplexity(), and setPreviousPeriodicComplexity().

◆ previousPosition_

Vec3D BaseParticle::previousPosition_
private

Indicates if this particle belongs to the maser boundary or is released into the wide open world.

Referenced by BaseParticle(), getPreviousPosition(), movePrevious(), and setPreviousPosition().

◆ radius_

◆ timeStamp_

unsigned BaseParticle::timeStamp_
private

Stores the time stamp the particle was created.

Referenced by BaseParticle(), getTimeStamp(), read(), setTimeStamp(), and write().


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