ParticleSpecies Class Referenceabstract

#include <ParticleSpecies.h>

+ Inheritance diagram for ParticleSpecies:

Public Types

typedef BaseInteraction InteractionType
 
typedef BaseSpecies MixedSpeciesType
 

Public Member Functions

 ParticleSpecies ()
 The default constructor. More...
 
 ParticleSpecies (const ParticleSpecies &p)
 The default copy constructor. More...
 
 ParticleSpecies (BaseNormalForce *normalForce, BaseFrictionForce *frictionForce, BaseAdhesiveForce *adhesiveForce)
 
 ~ParticleSpecies ()
 The default destructor. More...
 
ParticleSpeciescopy () const override=0
 Creates a deep copy of the object from which it is called. More...
 
virtual BaseSpeciescopyMixed () const =0
 Creates a new MixedSpecies with the same force properties as the Species from which it is called. See Species::copyMixed for details. More...
 
void read (std::istream &is) override
 Reads the species properties from an input stream. More...
 
void write (std::ostream &os) const override
 Writes the species properties to an output stream. More...
 
std::string getBaseName () const
 Used in Species::getName to obtain a unique name for each Species. More...
 
void setDensity (Mdouble density)
 
Mdouble getMassFromRadius (Mdouble radius) const
 
Mdouble getMassFromRadius (const Mdouble radius, SpeciesHandler &speciesHandler)
 
Mdouble getVolumeFromRadius (Mdouble radius) const
 
Mdouble getDensity () const
 Allows density_ to be accessed. More...
 
void computeMass (BaseParticle *p) const
 Compute Particle mass function, which required a reference to the Species vector. It computes the Particles mass, Inertia and the inverses. More...
 
void setTemperatureDependentDensity (const std::function< double(double)> &temperatureDependentDensity)
 
const std::function< double(double)> & getTemperatureDependentDensity () const
 
Mdouble getSmallestParticleMass () const
 Computes mass of the lightest particle (by mass) belonging to this species. This computation calls getLightestInverseParticleMassLocal, such that the computation is done on each node. More...
 
Mdouble getMaxInteractionDistance () const
 returns the largest separation distance at which adhesive short-range forces can occur. More...
 
void setMaxInteractionDistance (Mdouble interactionDistance=0)
 Sets maxInteractionDistance_. More...
 
const BaseSpeciesgetMixedSpecies (const ParticleSpecies *s) const
 
virtual void actionsAfterTimeStep (BaseParticle *particle) const
 
- Public Member Functions inherited from BaseSpecies
 BaseSpecies ()
 The default constructor. More...
 
 BaseSpecies (BaseNormalForce *normalForce, BaseFrictionForce *frictionForce_, BaseAdhesiveForce *adhesiveForce)
 
 BaseSpecies (const BaseSpecies &p)
 The copy constructor. More...
 
 ~BaseSpecies ()
 The default destructor. More...
 
virtual void copyInto (BaseSpecies *s) const =0
 
void setHandler (SpeciesHandler *handler)
 Sets the pointer to the handler to which this species belongs. More...
 
SpeciesHandlergetHandler () const
 Returns the pointer to the handler to which this species belongs. More...
 
virtual void mixAll (BaseSpecies *S, BaseSpecies *T)=0
 creates default values for mixed species More...
 
virtual bool getUseAngularDOFs () const =0
 Returns true if torques (i.e. angular degrees of freedom) have to be calculated. More...
 
virtual BaseInteractiongetNewInteraction (BaseInteractable *P, BaseInteractable *I, unsigned timeStamp) const =0
 returns new Interaction object. More...
 
virtual BaseInteractiongetEmptyInteraction () const =0
 
virtual void deleteEmptyInteraction (BaseInteraction *interaction) const =0
 
Mdouble getInteractionDistance () const
 returns the largest separation distance at which adhesive short-range forces can occur. More...
 
BaseNormalForcegetNormalForce () const
 
BaseFrictionForcegetFrictionForce () const
 
BaseAdhesiveForcegetAdhesiveForce () const
 
void setInteractionDistance (Mdouble interactionDistance)
 
- 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 std::string getName () const =0
 A purely virtual function. 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
 

Private Member Functions

Mdouble getLargestInverseParticleMassLocal () const
 Computes inverse mass of the lightest particle (by mass) belonging to this species. If MPI is used, this computation is done locally on each node. More...
 

Private Attributes

Mdouble density_
 The mass density. More...
 
std::function< double(double temperature)> temperatureDependentDensity_
 
Mdouble maxInteractionDistance_
 

Additional Inherited Members

- Static Public Member Functions inherited from BaseSpecies
static Mdouble average (Mdouble a, Mdouble b)
 Returns the harmonic mean of two variables. More...
 
static Mdouble averageInf (Mdouble a, Mdouble b)
 Returns the harmonic mean of two variables, returning inf if either is inf. More...
 
- Protected Attributes inherited from BaseSpecies
BaseNormalForcenormalForce_
 A pointer to the normal force parameters \detail This pointer is used by the Interaction's to get a pointer to the species The pointer is set in the constructors of SPecies and MixedSpecies. More...
 
BaseFrictionForcefrictionForce_
 A pointer to the friction force parameters \detail This pointer is used by the Interaction's to get a pointer to the species The pointer is set in the constructors of SPecies and MixedSpecies. More...
 
BaseAdhesiveForceadhesiveForce_
 A pointer to the adhesive force parameters \detail This pointer is used by the Interaction's to get a pointer to the species The pointer is set in the constructors of SPecies and MixedSpecies. More...
 

Member Typedef Documentation

◆ InteractionType

◆ MixedSpeciesType

Constructor & Destructor Documentation

◆ ParticleSpecies() [1/3]

ParticleSpecies::ParticleSpecies ( )

The default constructor.

38  : BaseSpecies(), density_(1.0)
39 {
40 #ifdef DEBUG_CONSTRUCTOR
41  std::cout<<"ParticleSpecies::ParticleSpecies() finished"<<std::endl;
42 #endif
43 }
BaseSpecies()
The default constructor.
Definition: BaseSpecies.cc:38
Mdouble density_
The mass density.
Definition: ParticleSpecies.h:136

◆ ParticleSpecies() [2/3]

ParticleSpecies::ParticleSpecies ( const ParticleSpecies p)

The default copy constructor.

Parameters
[in]pthe species that is copied
58  : BaseSpecies(p)
59 {
60  density_ = p.density_;
62 #ifdef DEBUG_CONSTRUCTOR
63  std::cout<<"ParticleSpecies::ParticleSpecies(const ParticleSpecies &p) finished"<<std::endl;
64 #endif
65 }
std::function< double(double temperature)> temperatureDependentDensity_
Definition: ParticleSpecies.h:142

References density_, and temperatureDependentDensity_.

◆ ParticleSpecies() [3/3]

ParticleSpecies::ParticleSpecies ( BaseNormalForce normalForce,
BaseFrictionForce frictionForce,
BaseAdhesiveForce adhesiveForce 
)
46  : BaseSpecies(normalForce,frictionForce, adhesiveForce), density_(1.0)
47 {
48 #ifdef DEBUG_CONSTRUCTOR
49  std::cout<<"ParticleSpecies::ParticleSpecies(n,f,a) finished"<<std::endl;
50 #endif
51 }

◆ ~ParticleSpecies()

ParticleSpecies::~ParticleSpecies ( )

The default destructor.

68 {
69 #ifdef DEBUG_DESTRUCTOR
70  std::cout<<"ParticleSpecies::~ParticleSpecies() finished"<<std::endl;
71 #endif
72 }

Member Function Documentation

◆ actionsAfterTimeStep()

◆ computeMass()

void ParticleSpecies::computeMass ( BaseParticle p) const

Compute Particle mass function, which required a reference to the Species vector. It computes the Particles mass, Inertia and the inverses.

Compute BaseParticle mass function, which required a reference to the Species vector. It computes the Particles mass, Inertia and the inverses. this function is called, if BaseParticleHandler::addObject, SpeciesHandler::addObject, ParticleSpecies::setDensity, BaseParticle::setRadius or DPMBase::setParticleDimensions is called

168 {
169  p->computeMass(*this);
170 }
virtual void computeMass(const ParticleSpecies &s)
Computes the particle's (inverse) mass and inertia.
Definition: BaseParticle.cc:876

References BaseParticle::computeMass().

Referenced by SuperQuadricParticle::setAxes(), SuperQuadricParticle::setExponents(), BaseParticle::setRadius(), WallParticleCollision::setupInitialConditions(), test1(), test2(), and BaseParticle::unfix().

◆ copy()

◆ copyMixed()

virtual BaseSpecies* ParticleSpecies::copyMixed ( ) const
pure virtual

◆ getBaseName()

std::string ParticleSpecies::getBaseName ( ) const

Used in Species::getName to obtain a unique name for each Species.

Returns
a string containing the name of the species (minus the word "Species")
101 {
102  return "Particle";
103 }

◆ getDensity()

◆ getLargestInverseParticleMassLocal()

Mdouble ParticleSpecies::getLargestInverseParticleMassLocal ( ) const
private

Computes inverse mass of the lightest particle (by mass) belonging to this species. If MPI is used, this computation is done locally on each node.

Returns
A pointer to the to the lightest BaseParticle (by mass) in this ParticleHandler.
189 {
190  Mdouble maxInvMass = 0;
191  logger.assert_debug(getHandler() != nullptr && getHandler()->getDPMBase() != nullptr,"speciesHandler must be set");
192  for (BaseParticle* const p : getHandler()->getDPMBase()->particleHandler)
193  {
194  if (p->getSpecies()==this && !(p->isFixed() || p->isMPIParticle() || p->isPeriodicGhostParticle()) && p->getInvMass() > maxInvMass)
195  {
196  maxInvMass = p->getInvMass();
197  }
198  }
199  return maxInvMass;
200 }
double Mdouble
Definition: GeneralDefine.h:34
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
Definition: BaseParticle.h:54
SpeciesHandler * getHandler() const
Returns the pointer to the handler to which this species belongs.
Definition: BaseSpecies.cc:99

References BaseSpecies::getHandler(), and logger.

Referenced by getSmallestParticleMass().

◆ getMassFromRadius() [1/2]

Mdouble ParticleSpecies::getMassFromRadius ( const Mdouble  radius,
SpeciesHandler speciesHandler 
)
129 {
130  setHandler(&speciesHandler);
131  return getMassFromRadius(radius);
132 }
void setHandler(SpeciesHandler *handler)
Sets the pointer to the handler to which this species belongs.
Definition: BaseSpecies.cc:91
Mdouble getMassFromRadius(Mdouble radius) const
Definition: ParticleSpecies.cc:123

References getMassFromRadius(), and BaseSpecies::setHandler().

◆ getMassFromRadius() [2/2]

Mdouble ParticleSpecies::getMassFromRadius ( Mdouble  radius) const
Todo:
TW: should getMassFromRadius be removed? IFCD: it is used in at least one driver (AxisymmetricHopper).
124 {
125  return getDensity() * getVolumeFromRadius(radius);
126 }
Mdouble getDensity() const
Allows density_ to be accessed.
Definition: ParticleSpecies.cc:118
Mdouble getVolumeFromRadius(Mdouble radius) const
Definition: ParticleSpecies.cc:135

References getDensity(), and getVolumeFromRadius().

Referenced by NautaMixer::addSpeciesAndSetTimeStepAndSaveCount(), BaseCluster::calculateTimeStep(), DropletBoundary::checkBoundaryAfterParticlesMove(), ClosedCSCWalls::ClosedCSCWalls(), CSCWalls::CSCWalls(), getMassFromRadius(), MarbleRun::getParticleMass(), main(), ChuteBottom::makeRoughBottom(), ParticleCreation::ParticleCreation(), regimeForceUnitTest::regimeForceUnitTest(), BaseCluster::setupInitialConditions(), SinterPair::SinterPair(), and viscoElasticUnitTest::viscoElasticUnitTest().

◆ getMaxInteractionDistance()

Mdouble ParticleSpecies::getMaxInteractionDistance ( ) const
inline

returns the largest separation distance at which adhesive short-range forces can occur.

returns the largest separation distance (negative overlap) at which (adhesive) short-range forces can occur (needed for contact detection). Defined in each of the AdhesiveForceSpecies It is defined as a virtual function here to allow the function to be called from a BaseSpecies pointer (which is the kind of pointer used for MixedSpecies).

113 {return maxInteractionDistance_;}
Mdouble maxInteractionDistance_
Definition: ParticleSpecies.h:147

References maxInteractionDistance_.

Referenced by BaseParticle::getMaxInteractionRadius().

◆ getMixedSpecies()

const BaseSpecies * ParticleSpecies::getMixedSpecies ( const ParticleSpecies s) const
238  {
239  return (getIndex()==s->getIndex())?this:getHandler()->getMixedObject(getIndex(),s->getIndex());
240 }
unsigned int getIndex() const
Returns the index of the object in the handler.
Definition: BaseObject.h:118
std::enable_if<!std::is_pointer< typename U::MixedSpeciesType >::value, typename U::MixedSpeciesType * >::type getMixedObject(const U *S, const U *T)
Definition: SpeciesHandler.h:74

References BaseSpecies::getHandler(), BaseObject::getIndex(), and SpeciesHandler::getMixedObject().

Referenced by BaseParticle::getInteractionDistance().

◆ getSmallestParticleMass()

Mdouble ParticleSpecies::getSmallestParticleMass ( ) const

Computes mass of the lightest particle (by mass) belonging to this species. This computation calls getLightestInverseParticleMassLocal, such that the computation is done on each node.

203 {
204 #ifdef MERCURYDPM_USE_MPI
205  Mdouble maxInvMass = 0;
207  //Obtain the global value
208  MPIContainer::Instance().allReduce(invMassLocal, maxInvMass, MPI_MAX);
209  //return value
210  return 1.0 / maxInvMass;
211 #else
212  return 1.0 / getLargestInverseParticleMassLocal();
213 #endif
214 
215 }
static MPIContainer & Instance()
fetch the instance to be used for communication
Definition: MpiContainer.h:134
Mdouble getLargestInverseParticleMassLocal() const
Computes inverse mass of the lightest particle (by mass) belonging to this species....
Definition: ParticleSpecies.cc:188

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

◆ getTemperatureDependentDensity()

const std::function< double(double)> & ParticleSpecies::getTemperatureDependentDensity ( ) const
173 {
175 }

References temperatureDependentDensity_.

◆ getVolumeFromRadius()

Mdouble ParticleSpecies::getVolumeFromRadius ( Mdouble  radius) const
Todo:
this should depend on the particle shape; thus, it should be a static function of BaseParticle
136 {
137  if (getHandler() == nullptr)
138  {
139  logger(ERROR,
140  "[Species::VolumeFromRadius()] No handler has been set, therefore, I can't figure out the dimensions.");
141  return 0;
142  }
143 
144  unsigned int particleDimensions = getHandler()->getDPMBase()->getParticleDimensions();
145  if (particleDimensions == 3)
146  {
147  return 4.0 / 3.0 * constants::pi * radius * radius * radius;
148  }
149  else if (particleDimensions == 2)
150  {
151  return constants::pi * radius * radius;
152  }
153  else if (particleDimensions == 1)
154  {
155  return 2.0 * radius;
156  }
157  else
158  {
159  logger(ERROR, "[Species::VolumeFromRadius()] the dimension of the particle is wrongly set to %",
160  particleDimensions);
161  return 0.0;
162  }
163 }
@ ERROR
DPMBase * getDPMBase()
Gets the problem that is solved using this handler.
Definition: BaseHandler.h:725
unsigned int getParticleDimensions() const
Returns the particle dimensionality.
Definition: DPMBase.cc:1467
const Mdouble pi
Definition: ExtendedMath.h:45

References ERROR, BaseHandler< T >::getDPMBase(), BaseSpecies::getHandler(), DPMBase::getParticleDimensions(), logger, and constants::pi.

Referenced by getMassFromRadius().

◆ read()

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

Reads the species properties from an input stream.

Parameters
[in]isinput stream (typically the restart file)

Reimplemented from BaseSpecies.

Reimplemented in Species< NormalForceSpecies, FrictionForceSpecies, AdhesiveForceSpecies >, Species< LinearViscoelasticNormalSpecies >, Species< LinearPlasticViscoelasticNormalSpecies, SlidingFrictionSpecies, IrreversibleAdhesiveSpecies >, and Species< LinearViscoelasticNormalSpecies, SlidingFrictionSpecies, IrreversibleAdhesiveSpecies >.

90 {
91  BaseObject::read(is);
92  std::string dummy;
93  is >> dummy >> density_;
95 }
virtual void read(std::istream &is)=0
Definition: BaseObject.cc:81
void read(std::istream &is) override
Definition: BaseSpecies.cc:140

References density_, BaseSpecies::read(), and BaseObject::read().

Referenced by Species< NormalForceSpecies, FrictionForceSpecies, AdhesiveForceSpecies >::read().

◆ setDensity()

void ParticleSpecies::setDensity ( Mdouble  density)

Allows density_ to be changed

Todo:
recalculate masses when setting dim_particle or rho
Parameters
[in]densitythe particle density
109 {
110  logger.assert_always(density > 0, "[ParticleSpecies::setDensity(%)] value has to be positive", density);
111  density_ = density;
113 }
ParticleHandler particleHandler
An object of the class ParticleHandler, contains the pointers to all the particles created.
Definition: DPMBase.h:1437
void computeAllMasses(unsigned int indSpecies)
Computes the mass for all BaseParticle of the given species in this ParticleHandler.
Definition: ParticleHandler.cc:1215

References ParticleHandler::computeAllMasses(), density_, BaseHandler< T >::getDPMBase(), BaseSpecies::getHandler(), BaseObject::getIndex(), logger, and DPMBase::particleHandler.

Referenced by NautaMixer::addSpeciesAndSetTimeStepAndSaveCount(), AngledPeriodicBoundarySecondUnitTest::AngledPeriodicBoundarySecondUnitTest(), AngledPeriodicBoundaryUnitTest::AngledPeriodicBoundaryUnitTest(), AngleOfRepose::AngleOfRepose(), BouncingSuperQuadric::BouncingSuperQuadric(), BoundariesSelfTest::BoundariesSelfTest(), Membrane::buildMesh(), ChutePeriodicDemo::ChutePeriodicDemo(), ChuteWithWedge::ChuteWithWedge(), ClosedCSCWalls::ClosedCSCWalls(), Contact::Contact(), ParameterStudy1DDemo::createSpecies(), ParameterStudy2DDemo::createSpecies(), ParameterStudy3DDemo::createSpecies(), CSCWalls::CSCWalls(), FluxAndPeriodicBoundarySelfTest::FluxAndPeriodicBoundarySelfTest(), FluxBoundarySelfTest::FluxBoundarySelfTest(), ForceLawsMPI2Test::ForceLawsMPI2Test(), MembraneDemo::initializeSpecies(), MembraneSelfTest::initializeSpecies(), InsertionBoundaryMPI2Test::InsertionBoundaryMPI2Test(), InsertionBoundarySelfTest::InsertionBoundarySelfTest(), LiquidMigrationMPI2Test::LiquidMigrationMPI2Test(), main(), MaserRepeatedOutInMPI2Test::MaserRepeatedOutInMPI2Test(), MercuryCGSelfTest::MercuryCGSelfTest(), MinimalExampleDrum::MinimalExampleDrum(), MultiParticlesInsertion::MultiParticlesInsertion(), ParticleCreation::ParticleCreation(), ParticleInclusion::ParticleInclusion(), protectiveWall::protectiveWall(), FileReader::read(), SpeciesHandler::readOldObject(), regimeForceUnitTest::regimeForceUnitTest(), FlowRule::setDensityVariation(), MercuryOS::setMaterialProperties(), MarbleRun::setParticleDensity(), GranularJet::setSilbert(), ChutePeriodic::setup(), NozzleDemo::setupInitialConditions(), Drum::setupInitialConditions(), CubicCell::setupInitialConditions(), FreeCooling2DinWallsDemo::setupInitialConditions(), FreeCoolingDemoProblem::setupInitialConditions(), ShiftingConstantMassFlowMaserBoundarySelfTest::setupInitialConditions(), ShiftingMaserBoundarySelfTest::setupInitialConditions(), GetDistanceAndNormalForIntersectionOfWalls::setupInitialConditions(), GetDistanceAndNormalForScrew::setupInitialConditions(), GetDistanceAndNormalForTriangleWall::setupInitialConditions(), ConstantMassFlowMaserBoundaryMixedSpeciesSelfTest::setupInitialConditions(), ConstantMassFlowMaserSelfTest::setupInitialConditions(), CubeDeletionBoundarySelfTest::setupInitialConditions(), DeletionBoundarySelfTest::setupInitialConditions(), DistributionSelfTest::setupInitialConditions(), DistributionToPSDSelfTest::setupInitialConditions(), InsertionBoundarySelfTest::setupInitialConditions(), MultiplePSDSelfTest::setupInitialConditions(), NozzleSelfTest::setupInitialConditions(), PolydisperseInsertionBoundarySelfTest::setupInitialConditions(), PSDManualInsertionSelfTest::setupInitialConditions(), PSDSelfTest::setupInitialConditions(), StressStrainControl::setupInitialConditions(), SubcriticalMaserBoundarySelfTest::setupInitialConditions(), ParticleParticleCollision::setupInitialConditions(), WallParticleCollision::setupInitialConditions(), ContactDetectionIntersectionOfWallsTest::setupInitialConditions(), GetDistanceAndNormalForTriangleWalls::setupInitialConditions(), RollingOverTriangleWalls::setupInitialConditions(), UnionOfWalls::setupInitialConditions(), EllipsoidsBouncingOnWallDemo::setupInitialConditions(), EllipticalSuperQuadricCollision::setupInitialConditions(), SlidingSpheresUnitTest::setupInitialConditions(), ShapesDemo::setupInitialConditions(), ParticleParticleInteraction::setupInitialConditions(), ParticleParticleInteractionWithPlasticForces::setupInitialConditions(), ParticleWallInteraction::setupInitialConditions(), Packing::setupInitialConditions(), CreateDataAndFStatFiles::setupInitialConditions(), RandomClusterInsertionBoundarySelfTest::setupInitialConditions(), FreeFall::setupInitialConditions(), HertzianSinterForceUnitTest::setupInitialConditions(), MovingIntersectionOfWallsUnitTest_Basic::setupInitialConditions(), MovingWalls::setupInitialConditions(), MovingWall::setupInitialConditions(), PeriodicWallsWithSlidingFrictionUnitTest::setupInitialConditions(), PlasticForceUnitTest::setupInitialConditions(), SeparateFilesSelfTest::setupInitialConditions(), WallSpecies::setupInitialConditions(), Siegen::Siegen(), SilbertPeriodic::SilbertPeriodic(), SinterPair::SinterPair(), StressStrainControl::StressStrainControl(), SubcriticalMaserBoundaryTESTMPI2Test::SubcriticalMaserBoundaryTESTMPI2Test(), T_protectiveWall::T_protectiveWall(), test1(), test2(), viscoElasticUnitTest::viscoElasticUnitTest(), and Wall::Wall().

◆ setMaxInteractionDistance()

void ParticleSpecies::setMaxInteractionDistance ( Mdouble  interactionDistance = 0)

Sets maxInteractionDistance_.

Parameters
interactionDistancethe interaction distance that has been changed

Sets maxInteractionDistance_

Parameters
interactionDistancethe interaction distance that has been changed
221  {
222  // if maxInteractionDistance_ has increased it's simple
223  if (interactionDistance>=maxInteractionDistance_) {
224  maxInteractionDistance_ = interactionDistance;
225  } else /*else we need to recompute*/ {
227  int j = getIndex();
228  for (int i=0; i<getHandler()->getSize(); ++i) {
229  const auto mixedSpecies = getHandler()->getMixedObject(i,j);
230  if (mixedSpecies) { //this check is necessary because the mixed species handler might not yet be built fully
232  mixedSpecies->getInteractionDistance());
233  }
234  }
235  }
236 }
unsigned int getSize() const
Gets the size of the particleHandler (including mpi and periodic particles)
Definition: BaseHandler.h:655
const std::complex< Mdouble > i
Definition: ExtendedMath.h:51

References BaseSpecies::getHandler(), BaseObject::getIndex(), SpeciesHandler::getMixedObject(), BaseHandler< T >::getSize(), constants::i, and maxInteractionDistance_.

Referenced by SpeciesHandler::addObject(), and BaseSpecies::setInteractionDistance().

◆ setTemperatureDependentDensity()

void ParticleSpecies::setTemperatureDependentDensity ( const std::function< double(double)> &  temperatureDependentDensity)
179 {
180  temperatureDependentDensity_ = temperatureDependentDensity;
181 // density_ = temperatureDependentDensity_(0);
182 // logger(INFO,"Setting initial density to %",temperature_);
183 }

References temperatureDependentDensity_.

◆ write()

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

Writes the species properties to an output stream.

Parameters
[out]osoutput stream (typically the restart file)

Reimplemented from BaseSpecies.

Reimplemented in Species< NormalForceSpecies, FrictionForceSpecies, AdhesiveForceSpecies >, Species< LinearViscoelasticNormalSpecies >, Species< LinearPlasticViscoelasticNormalSpecies, SlidingFrictionSpecies, IrreversibleAdhesiveSpecies >, and Species< LinearViscoelasticNormalSpecies, SlidingFrictionSpecies, IrreversibleAdhesiveSpecies >.

78 {
79  //note we inherit from BaseObject, not BaseParticle
81  os << " density " << density_;
83  //todo flip the two values
84 }
virtual void write(std::ostream &os) const =0
A purely virtual function which has an implementation which writes the name and the object id_ to the...
Definition: BaseObject.cc:91
void write(std::ostream &os) const override
Sets the boolean constantRestitution_.
Definition: BaseSpecies.cc:131

References density_, BaseObject::write(), and BaseSpecies::write().

Referenced by Species< NormalForceSpecies, FrictionForceSpecies, AdhesiveForceSpecies >::write().

Member Data Documentation

◆ density_

Mdouble ParticleSpecies::density_
private

The mass density.

Referenced by getDensity(), ParticleSpecies(), read(), setDensity(), and write().

◆ maxInteractionDistance_

Mdouble ParticleSpecies::maxInteractionDistance_
private

Returns the max distance between particles of this species and any other species below which adhesive forces can occur (needed for contact detection)

Referenced by getMaxInteractionDistance(), and setMaxInteractionDistance().

◆ temperatureDependentDensity_

std::function<double(double temperature)> ParticleSpecies::temperatureDependentDensity_
private

Change this function to let the particles expand due to temperature. The default value (empty) stands for constant density.

Referenced by getTemperatureDependentDensity(), ParticleSpecies(), and setTemperatureDependentDensity().


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