SPHNormalSpecies Class Reference

SPHNormalSpecies contains the parameters used to describe a SPH model. More...

#include <SPHNormalSpecies.h>

+ Inheritance diagram for SPHNormalSpecies:

Public Types

typedef SPHInteraction InteractionType
 The correct Interaction type for this FrictionForceSpecies. More...
 

Public Member Functions

 SPHNormalSpecies ()
 The default constructor. More...
 
 SPHNormalSpecies (const SPHNormalSpecies &p)
 The default copy constructor. More...
 
 ~SPHNormalSpecies ()
 The default destructor. More...
 
void read (std::istream &is)
 Reads the species properties from an input stream. More...
 
void write (std::ostream &os) const
 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...
 
Mdouble getMaximumVelocity (Mdouble radius, Mdouble mass) const
 Calculates the maximum velocity allowed for a collision of two copies of P (for higher velocities particles could pass through each other) More...
 
void setStiffnessAndRestitutionCoefficient (Mdouble k_, Mdouble eps, Mdouble mass)
 Sets k, disp such that it matches a given tc and eps for a collision of two copies of P. More...
 
void setRestitutionCoefficient (double eps, Mdouble mass)
 Sets disp to obtain a restitution coefficient eps for a collision of two particles of mass m. More...
 
void setCollisionTimeAndRestitutionCoefficient (Mdouble tc, Mdouble eps, BaseParticle *p)
 Sets k, disp such that it matches a given tc and eps for a collision of two copies of particle p. More...
 
void setCollisionTimeAndRestitutionCoefficient (Mdouble tc, Mdouble eps, Mdouble mass)
 Sets k, disp such that it matches a given tc and eps for a collision of two copies of equal mass m. More...
 
void setCollisionTimeAndRestitutionCoefficient (Mdouble collisionTime, Mdouble restitutionCoefficient, Mdouble mass1, Mdouble mass2)
 
Mdouble getCollisionTime (Mdouble mass) const
 Calculates collision time for two copies of a particle of given disp, k, mass. More...
 
Mdouble getRestitutionCoefficient (Mdouble mass) const
 Calculates restitution coefficient for two copies of given disp, k, mass. More...
 
void mix (SPHNormalSpecies *SBase, SPHNormalSpecies *TBase)
 creates default values for mixed species More...
 
void setStiffness (Mdouble new_k)
 Allows the spring constant to be changed. More...
 
Mdouble getStiffness () const
 Allows the spring constant to be accessed. More...
 
void setDissipation (Mdouble dissipation)
 Allows the normal dissipation to be changed. More...
 
Mdouble getDissipation () const
 Allows the normal dissipation to be accessed. More...
 
MERCURYDPM_DEPRECATED void setStiffnessAndDissipation (helpers::KAndDisp new_)
 Allows the spring and dissipation constants to be changed simultaneously. More...
 
- Public Member Functions inherited from BaseNormalForce
 BaseNormalForce ()
 
 BaseNormalForce (const BaseNormalForce &p)
 
bool getConstantRestitution () const
 
void setConstantRestitution (bool constantRestitution)
 
virtual void actionsAfterTimeStep (BaseParticle *particle) const
 
- Public Member Functions inherited from BaseForce
BaseSpeciesgetBaseSpecies () const
 
void setBaseSpecies (BaseSpecies *baseSpecies)
 

Private Attributes

Mdouble stiffness_
 (normal) spring constant More...
 
Mdouble dissipation_
 (normal) viscosity More...
 

Detailed Description

SPHNormalSpecies contains the parameters used to describe a SPH model.

See SPHNormalInteraction::computeForce for a description of the force law.

Member Typedef Documentation

◆ InteractionType

The correct Interaction type for this FrictionForceSpecies.

Constructor & Destructor Documentation

◆ SPHNormalSpecies() [1/2]

SPHNormalSpecies::SPHNormalSpecies ( )

The default constructor.

36  : BaseNormalForce()
37 {
38  stiffness_ = 0.0;
39  dissipation_ = 0.0;
40 #ifdef DEBUG_CONSTRUCTOR
41  std::cout<<"LinearViscoelasticNormalSpecies::LinearViscoelasticNormalSpecies() finished"<<std::endl;
42 #endif
43 }
BaseNormalForce()
Definition: BaseNormalForce.h:35
Mdouble dissipation_
(normal) viscosity
Definition: SPHNormalSpecies.h:117
Mdouble stiffness_
(normal) spring constant
Definition: SPHNormalSpecies.h:116

References dissipation_, and stiffness_.

◆ SPHNormalSpecies() [2/2]

SPHNormalSpecies::SPHNormalSpecies ( const SPHNormalSpecies p)

The default copy constructor.

Parameters
[in]thespecies that is copied
49  : BaseNormalForce(p)
50 {
53 #ifdef DEBUG_CONSTRUCTOR
54  std::cout<<"LinearViscoelasticNormalSpecies::LinearViscoelasticNormalSpecies(const LinearViscoelasticNormalSpecies &p) finished"<<std::endl;
55 #endif
56 }

References dissipation_, and stiffness_.

◆ ~SPHNormalSpecies()

SPHNormalSpecies::~SPHNormalSpecies ( )

The default destructor.

59 {
60 #ifdef DEBUG_DESTRUCTOR
61  std::cout<<"LinearViscoelasticNormalSpecies::~LinearViscoelasticNormalSpecies() finished"<<std::endl;
62 #endif
63 }

Member Function Documentation

◆ getBaseName()

std::string SPHNormalSpecies::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")
88 {
89  return "LinearViscoelastic";
90 }

◆ getCollisionTime()

Mdouble SPHNormalSpecies::getCollisionTime ( Mdouble  mass) const

Calculates collision time for two copies of a particle of given disp, k, mass.

Calculates collision time for two copies of a particle of given disp, k, mass

Parameters
[in]massmass of a typical particle
140 {
141  if (mass <= 0)
142  {
143  std::cerr << "Warning in getCollisionTime(" << mass
144  << ") mass is not set or has an unexpected value, (getCollisionTime(" << mass << "))" << std::endl;
145  }
146  if (stiffness_ <= 0)
147  {
148  std::cerr << "Warning in getCollisionTime(" << mass << ") stiffness=" << stiffness_
149  << " is not set or has an unexpected value, (getCollisionTime(" << mass << "), with stiffness="
150  << stiffness_ << ")" << std::endl;
151  }
152  if (dissipation_ < 0)
153  {
154  std::cerr << "Warning in getCollisionTime(" << mass << ") dissipation=" << dissipation_
155  << " is not set or has an unexpected value, (getCollisionTime(" << mass << "), with dissipation="
156  << dissipation_ << ")" << std::endl;
157  }
158  Mdouble tosqrt = stiffness_ / (.5 * mass) - mathsFunc::square(dissipation_ / mass);
159  if (tosqrt <= 0)
160  {
161  std::cerr << "Warning in getCollisionTime(" << mass
162  << ") values for mass, stiffness and dissipation would lead to an overdamped system, (getCollisionTime("
163  << mass << "), with stiffness=" << stiffness_ << " and dissipation=" << dissipation_ << ")"
164  << std::endl;
165  }
166  return constants::pi / std::sqrt(tosqrt);
167 }
double Mdouble
Definition: GeneralDefine.h:34
const Mdouble pi
Definition: ExtendedMath.h:45
T square(const T val)
squares a number
Definition: ExtendedMath.h:106

References dissipation_, constants::pi, mathsFunc::square(), and stiffness_.

Referenced by getRestitutionCoefficient().

◆ getDissipation()

Mdouble SPHNormalSpecies::getDissipation ( ) const

Allows the normal dissipation to be accessed.

133 {
134  return dissipation_;
135 }

References dissipation_.

Referenced by SPHInteraction::computeNormalForce(), mix(), and setCollisionTimeAndRestitutionCoefficient().

◆ getMaximumVelocity()

Mdouble SPHNormalSpecies::getMaximumVelocity ( Mdouble  radius,
Mdouble  mass 
) const

Calculates the maximum velocity allowed for a collision of two copies of P (for higher velocities particles could pass through each other)

177 {
178  return radius * std::sqrt(stiffness_ / (.5 * mass));
179 }

References stiffness_.

◆ getRestitutionCoefficient()

Mdouble SPHNormalSpecies::getRestitutionCoefficient ( Mdouble  mass) const

Calculates restitution coefficient for two copies of given disp, k, mass.

171 {
172  return std::exp(-dissipation_ / mass * getCollisionTime(mass));
173 }
Mdouble getCollisionTime(Mdouble mass) const
Calculates collision time for two copies of a particle of given disp, k, mass.
Definition: SPHNormalSpecies.cc:139
Mdouble exp(Mdouble Exponent)
Definition: ExtendedMath.cc:84

References dissipation_, mathsFunc::exp(), and getCollisionTime().

◆ getStiffness()

Mdouble SPHNormalSpecies::getStiffness ( ) const

◆ mix()

void SPHNormalSpecies::mix ( SPHNormalSpecies SBase,
SPHNormalSpecies TBase 
)

creates default values for mixed species

For all parameters we assume that the harmonic mean of the parameters of the original two species is a sensible default.

Parameters
[in]S,Tthe two species whose properties are mixed to create the new species
263 {
264  stiffness_ = BaseSpecies::average(S->getStiffness(), T->getStiffness());
265  dissipation_ = BaseSpecies::average(S->getDissipation(), T->getDissipation());
266 }
static Mdouble average(Mdouble a, Mdouble b)
Returns the harmonic mean of two variables.
Definition: BaseSpecies.cc:110

References BaseSpecies::average(), dissipation_, getDissipation(), getStiffness(), and stiffness_.

◆ read()

void SPHNormalSpecies::read ( std::istream &  is)

Reads the species properties from an input stream.

Parameters
[in]inputstream (typically the restart file)
78 {
79  std::string dummy;
80  is >> dummy >> stiffness_
81  >> dummy >> dissipation_;
82 }

References dissipation_, and stiffness_.

◆ setCollisionTimeAndRestitutionCoefficient() [1/3]

void SPHNormalSpecies::setCollisionTimeAndRestitutionCoefficient ( Mdouble  collisionTime,
Mdouble  restitutionCoefficient,
Mdouble  mass1,
Mdouble  mass2 
)

Set k, disp such that is matches a given tc and eps for a collision of two different masses. Recall the resitution constant is a function of k, disp and the mass of each particle in the collision See also setCollisionTimeAndRestitutionCoefficient(Mdouble tc, Mdouble eps, Mdouble mass)

Set k, disp such that is matches a given tc and eps for a collision of two different masses. Recall the resitution constant is a function of k, disp and the mass of each particle in the collision See also setCollisionTimeAndRestitutionCoefficient(Mdouble tc, Mdouble eps, Mdouble mass)

Parameters
[in]collisiontime
251 {
252  Mdouble reduced_mass = mass1 * mass2 / (mass1 + mass2);
253  setCollisionTimeAndRestitutionCoefficient(collisionTime, restitutionCoefficient, 2.0 * reduced_mass);
254 }
void setCollisionTimeAndRestitutionCoefficient(Mdouble tc, Mdouble eps, BaseParticle *p)
Sets k, disp such that it matches a given tc and eps for a collision of two copies of particle p.
Definition: SPHNormalSpecies.cc:211

References setCollisionTimeAndRestitutionCoefficient().

◆ setCollisionTimeAndRestitutionCoefficient() [2/3]

void SPHNormalSpecies::setCollisionTimeAndRestitutionCoefficient ( Mdouble  tc,
Mdouble  eps,
BaseParticle p 
)

Sets k, disp such that it matches a given tc and eps for a collision of two copies of particle p.

212 {
213  Mdouble mass = p->getSpecies()->getDensity() * p->getVolume();
215 }
const ParticleSpecies * getSpecies() const
Returns a pointer to the species of this BaseInteractable.
Definition: BaseInteractable.h:108
virtual Mdouble getVolume() const
Get Particle volume function, which required a reference to the Species vector. It returns the volume...
Definition: BaseParticle.cc:143
Mdouble getDensity() const
Allows density_ to be accessed.
Definition: ParticleSpecies.cc:118

References ParticleSpecies::getDensity(), BaseInteractable::getSpecies(), and BaseParticle::getVolume().

Referenced by setCollisionTimeAndRestitutionCoefficient().

◆ setCollisionTimeAndRestitutionCoefficient() [3/3]

void SPHNormalSpecies::setCollisionTimeAndRestitutionCoefficient ( Mdouble  tc,
Mdouble  eps,
Mdouble  mass 
)

Sets k, disp such that it matches a given tc and eps for a collision of two copies of equal mass m.

Sets k, disp such that it matches a given tc and eps for a collision of two copies of equal mass m

Parameters
[in]tccollision time
[in]epsrestitution coefficient
[in]massharmonic average particle mass, \(\frac{2}{1/m1+1/m2}\)
225 {
226  if (eps == 0.0)
227  {
228  stiffness_ = .5 * mass * mathsFunc::square(constants::pi / tc);
229  dissipation_ = std::sqrt(2.0 * mass * stiffness_);
230  }
231  else
232  {
233  dissipation_ = -mass / tc * mathsFunc::log(eps);
234  stiffness_ = .5 * mass * (mathsFunc::square(constants::pi / tc)
235  + mathsFunc::square(dissipation_ / mass));
236  }
237  logger(INFO,
238  "setCollisionTimeAndRestitutionCoefficient: set stiffness to % and dissipation to % to obtain a collision time of % and a restitution coefficient of % for two particles of mass %",
239  getStiffness(), getDissipation(), tc, eps, mass);
240 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ INFO
Mdouble getStiffness() const
Allows the spring constant to be accessed.
Definition: SPHNormalSpecies.cc:105
Mdouble getDissipation() const
Allows the normal dissipation to be accessed.
Definition: SPHNormalSpecies.cc:132
Mdouble log(Mdouble Power)
Definition: ExtendedMath.cc:104

References dissipation_, getDissipation(), getStiffness(), INFO, mathsFunc::log(), logger, constants::pi, mathsFunc::square(), and stiffness_.

◆ setDissipation()

void SPHNormalSpecies::setDissipation ( Mdouble  dissipation)

Allows the normal dissipation to be changed.

119 {
120  if (dissipation >= 0)
121  {
122  dissipation_ = dissipation;
123  }
124  else
125  {
126  std::cerr << "Error in setDissipation(" << dissipation << ")" << std::endl;
127  exit(-1);
128  }
129 }

References dissipation_.

Referenced by setStiffnessAndDissipation().

◆ setRestitutionCoefficient()

void SPHNormalSpecies::setRestitutionCoefficient ( double  eps,
Mdouble  mass 
)

Sets disp to obtain a restitution coefficient eps for a collision of two particles of mass m.

Sets k, disp such that it matches a given tc and eps for a collision of two copies of P

Parameters
[in]stiffnessstiffness
[in]epsrestitution coefficient
[in]masseffective particle mass, \(\frac{2}{1/m1+1/m2}\)
200 {
201  if (eps == 0.0) {
202  dissipation_ = std::sqrt(2.0 * mass * getStiffness());
203  } else {
204  const Mdouble logEps = log(eps);
205  dissipation_ = -std::sqrt(2.0 * mass * getStiffness()
206  / (constants::sqr_pi + mathsFunc::square(logEps))) * logEps;
207  }
208 }
const Mdouble sqr_pi
Definition: ExtendedMath.h:47

References dissipation_, getStiffness(), mathsFunc::log(), constants::sqr_pi, and mathsFunc::square().

Referenced by setStiffnessAndRestitutionCoefficient().

◆ setStiffness()

void SPHNormalSpecies::setStiffness ( Mdouble  new_k)

Allows the spring constant to be changed.

94 {
95  if (new_k >= 0)
96  stiffness_ = new_k;
97  else
98  {
99  std::cerr << "Error in set_k" << std::endl;
100  exit(-1);
101  }
102 }

References stiffness_.

Referenced by setStiffnessAndDissipation().

◆ setStiffnessAndDissipation()

void SPHNormalSpecies::setStiffnessAndDissipation ( helpers::KAndDisp  new_)

Allows the spring and dissipation constants to be changed simultaneously.

112 {
113  setStiffness(new_.k);
114  setDissipation(new_.disp);
115 }
void setDissipation(Mdouble dissipation)
Allows the normal dissipation to be changed.
Definition: SPHNormalSpecies.cc:118
void setStiffness(Mdouble new_k)
Allows the spring constant to be changed.
Definition: SPHNormalSpecies.cc:93
Mdouble k
Definition: FormulaHelpers.h:39
Mdouble disp
Definition: FormulaHelpers.h:40

References helpers::KAndDisp::disp, helpers::KAndDisp::k, setDissipation(), and setStiffness().

◆ setStiffnessAndRestitutionCoefficient()

void SPHNormalSpecies::setStiffnessAndRestitutionCoefficient ( Mdouble  stiffness,
Mdouble  eps,
Mdouble  mass 
)

Sets k, disp such that it matches a given tc and eps for a collision of two copies of P.

Sets k, disp such that it matches a given tc and eps for a collision of two copies of P

Parameters
[in]stiffnessstiffness
[in]epsrestitution coefficient
[in]harmonicmean of particle masses, \(\frac{2}{1/m1+1/m2}\)
188 {
189  stiffness_ = stiffness;
190  setRestitutionCoefficient(eps, mass);
191 }
void setRestitutionCoefficient(double eps, Mdouble mass)
Sets disp to obtain a restitution coefficient eps for a collision of two particles of mass m.
Definition: SPHNormalSpecies.cc:199

References setRestitutionCoefficient(), and stiffness_.

◆ write()

void SPHNormalSpecies::write ( std::ostream &  os) const

Writes the species properties to an output stream.

Parameters
[out]outputstream (typically the restart file)
69 {
70  os << " stiffness " << stiffness_
71  << " dissipation " << dissipation_;
72 }

References dissipation_, and stiffness_.

Member Data Documentation

◆ dissipation_

◆ stiffness_


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