LinearPlasticViscoelasticNormalSpecies Class Reference

LinearPlasticViscoelasticNormalSpecies contains the parameters used to describe a plastic-cohesive normal force (Stefan Ludings plastic-cohesive force model). More...

#include <LinearPlasticViscoelasticNormalSpecies.h>

+ Inheritance diagram for LinearPlasticViscoelasticNormalSpecies:

Public Types

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

Public Member Functions

 LinearPlasticViscoelasticNormalSpecies ()
 The default constructor. More...
 
 LinearPlasticViscoelasticNormalSpecies (const LinearPlasticViscoelasticNormalSpecies &p)
 The default copy constructor. More...
 
 ~LinearPlasticViscoelasticNormalSpecies ()
 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...
 
void mix (LinearPlasticViscoelasticNormalSpecies *S, LinearPlasticViscoelasticNormalSpecies *T)
 creates default values for mixed species More...
 
void setCollisionTimeAndRestitutionCoefficient (Mdouble tc, Mdouble eps, Mdouble mass)
 
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...
 
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 setPlasticParameters (Mdouble loadingStiffness, Mdouble unloadingStiffnessMax, Mdouble cohesionStiffness, Mdouble penetrationDepthMax)
 Sets all parameters of the linear plastic-viscoelastic normal force at once. More...
 
Mdouble getLoadingStiffness () const
 Returns the loading stiffness of the linear plastic-viscoelastic normal force. More...
 
Mdouble getUnloadingStiffnessMax () const
 Returns the maximum unloading stiffness of the linear plastic-viscoelastic normal force. More...
 
Mdouble getCohesionStiffness () const
 Returns the cohesive stiffness of the linear plastic-viscoelastic normal force. More...
 
Mdouble getPenetrationDepthMax () const
 Returns the maximum penetration depth of the linear plastic-viscoelastic normal force. More...
 
void setLoadingStiffness (Mdouble loadingStiffness)
 Sets the loading stiffness of the linear plastic-viscoelastic normal force. More...
 
void setUnloadingStiffnessMax (Mdouble unloadingStiffnessMax)
 Sets the maximum unloading stiffness of the linear plastic-viscoelastic normal force. More...
 
void setCohesionStiffness (Mdouble cohesionStiffness)
 Sets the cohesive stiffness of the linear plastic-viscoelastic normal force. More...
 
void setPenetrationDepthMax (Mdouble penetrationDepthMax)
 Sets the maximum penetration depth of the linear plastic-viscoelastic normal force. More...
 
void setDissipation (Mdouble dissipation)
 Sets the linear dissipation coefficient of the linear plastic-viscoelastic normal force. More...
 
MERCURYDPM_DEPRECATED void setLoadingStiffnessAndDissipation (helpers::KAndDisp new_)
 Allows the spring and dissipation constants to be changed simultaneously. More...
 
Mdouble computeTimeStep (Mdouble mass)
 Returns the optimal time step to resolve a collision of two particles of a given mass. More...
 
Mdouble getDissipation () const
 Allows the normal dissipation to be accessed. More...
 
Mdouble computeBondNumberMax (Mdouble harmonicMeanRadius, Mdouble gravitationalAcceleration) const
 1) Computes the maximum plastic overlap delta_p* = phi*r 2) Computes the overlap at which the maximum adhesive force is generated: delta_c* = delta_p* / (1+k_c/k_2*) 3) Computes the maximum adhesive force f_c* = k_c * delta_c* 4) Computes the maximum bond number Bo* = f_c* / (m*g) More...
 
bool getDoConstantUnloadingStiffness () const
 
void setDoConstantUnloadingStiffness (bool doConstantUnloadingStiffness)
 
- 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 loadingStiffness_
 (normal) spring constant (k_1) More...
 
Mdouble unloadingStiffnessMax_
 the maximum elastic constant (k_2^max) for plastic deformations More...
 
Mdouble cohesionStiffness_
 the adhesive spring constant (k^c) for plastic deformations More...
 
Mdouble penetrationDepthMax_
 the depth (relative to the normalized radius) at which k_2^max is used (phi_f) More...
 
Mdouble dissipation_
 linear dissipation coefficient More...
 
bool doConstantUnloadingStiffness_ = false
 

Detailed Description

LinearPlasticViscoelasticNormalSpecies contains the parameters used to describe a plastic-cohesive normal force (Stefan Ludings plastic-cohesive force model).

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

Member Typedef Documentation

◆ InteractionType

The correct Interaction type for this FrictionForceSpecies.

Constructor & Destructor Documentation

◆ LinearPlasticViscoelasticNormalSpecies() [1/2]

LinearPlasticViscoelasticNormalSpecies::LinearPlasticViscoelasticNormalSpecies ( )

The default constructor.

37  : BaseNormalForce()
38 {
39  loadingStiffness_ = 0.0;
41  cohesionStiffness_ = 0.0;
43  dissipation_ = 0.0;
45 #ifdef DEBUG_CONSTRUCTOR
46  std::cout<<"LinearPlasticViscoelasticNormalSpecies::LinearPlasticViscoelasticNormalSpecies() finished"<<std::endl;
47 #endif
48 }
BaseNormalForce()
Definition: BaseNormalForce.h:35
Mdouble loadingStiffness_
(normal) spring constant (k_1)
Definition: LinearPlasticViscoelasticNormalSpecies.h:175
Mdouble cohesionStiffness_
the adhesive spring constant (k^c) for plastic deformations
Definition: LinearPlasticViscoelasticNormalSpecies.h:181
Mdouble unloadingStiffnessMax_
the maximum elastic constant (k_2^max) for plastic deformations
Definition: LinearPlasticViscoelasticNormalSpecies.h:178
Mdouble dissipation_
linear dissipation coefficient
Definition: LinearPlasticViscoelasticNormalSpecies.h:187
bool doConstantUnloadingStiffness_
Definition: LinearPlasticViscoelasticNormalSpecies.h:190
Mdouble penetrationDepthMax_
the depth (relative to the normalized radius) at which k_2^max is used (phi_f)
Definition: LinearPlasticViscoelasticNormalSpecies.h:184

References cohesionStiffness_, dissipation_, doConstantUnloadingStiffness_, loadingStiffness_, penetrationDepthMax_, and unloadingStiffnessMax_.

◆ LinearPlasticViscoelasticNormalSpecies() [2/2]

LinearPlasticViscoelasticNormalSpecies::LinearPlasticViscoelasticNormalSpecies ( const LinearPlasticViscoelasticNormalSpecies p)

The default copy constructor.

Parameters
[in]thespecies that is copied
55  : BaseNormalForce(p)
56 {
63 #ifdef DEBUG_CONSTRUCTOR
64  std::cout<<"LinearPlasticViscoelasticNormalSpecies::LinearPlasticViscoelasticNormalSpecies(const LinearPlasticViscoelasticNormalSpecies &p) finished"<<std::endl;
65 #endif
66 }

References cohesionStiffness_, dissipation_, doConstantUnloadingStiffness_, loadingStiffness_, penetrationDepthMax_, and unloadingStiffnessMax_.

◆ ~LinearPlasticViscoelasticNormalSpecies()

LinearPlasticViscoelasticNormalSpecies::~LinearPlasticViscoelasticNormalSpecies ( )

The default destructor.

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

Member Function Documentation

◆ computeBondNumberMax()

Mdouble LinearPlasticViscoelasticNormalSpecies::computeBondNumberMax ( Mdouble  harmonicMeanRadius,
Mdouble  gravitationalAcceleration 
) const

1) Computes the maximum plastic overlap delta_p* = phi*r 2) Computes the overlap at which the maximum adhesive force is generated: delta_c* = delta_p* / (1+k_c/k_2*) 3) Computes the maximum adhesive force f_c* = k_c * delta_c* 4) Computes the maximum bond number Bo* = f_c* / (m*g)

339  {
340  if (getConstantRestitution()) {
341  //harmonicMeanRadius unused
342  const Mdouble plasticOverlapMax = penetrationDepthMax_;
343  const Mdouble overlapMaxCohesion = plasticOverlapMax / (1 + cohesionStiffness_ / unloadingStiffnessMax_);
344  const Mdouble cohesionAccelerationMax = cohesionStiffness_ * overlapMaxCohesion;
345  return cohesionAccelerationMax / gravitationalAcceleration;
346  } else {
347  const Mdouble plasticOverlapMax = penetrationDepthMax_ * harmonicMeanRadius;
348  const Mdouble overlapMaxCohesion = plasticOverlapMax / (1 + cohesionStiffness_ / unloadingStiffnessMax_);
349  const Mdouble cohesionForceMax = cohesionStiffness_ * overlapMaxCohesion;
350  auto species = dynamic_cast<const ParticleSpecies*>(getBaseSpecies());
351  logger.assert_debug(species,"computeBondNumberMax: species needs to be a ParticleSpecies");
352  const Mdouble gravitationalForce = gravitationalAcceleration * species->getMassFromRadius(harmonicMeanRadius);
353  return cohesionForceMax / gravitationalForce;
354  }
355 }
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.
BaseSpecies * getBaseSpecies() const
Definition: BaseForce.h:35
bool getConstantRestitution() const
Definition: BaseNormalForce.h:46
Definition: ParticleSpecies.h:37

References cohesionStiffness_, BaseForce::getBaseSpecies(), BaseNormalForce::getConstantRestitution(), logger, penetrationDepthMax_, and unloadingStiffnessMax_.

◆ computeTimeStep()

Mdouble LinearPlasticViscoelasticNormalSpecies::computeTimeStep ( Mdouble  mass)

Returns the optimal time step to resolve a collision of two particles of a given mass.

Calculates collision time for stiffest spring constant, divides by 50

Parameters
[in]theoptimal time step is computed to resolve a collision of two particles of this mass. If constant restitution is enabled, the collision time is mass-independent.
220 {
221  if (getConstantRestitution()) mass = 1;
222  return 0.02 * constants::pi /
223  std::sqrt(unloadingStiffnessMax_ / (.5 * mass) - mathsFunc::square(dissipation_ / mass));
224 }
const Mdouble pi
Definition: ExtendedMath.h:45
T square(const T val)
squares a number
Definition: ExtendedMath.h:106

References dissipation_, BaseNormalForce::getConstantRestitution(), constants::pi, mathsFunc::square(), and unloadingStiffnessMax_.

◆ getBaseName()

std::string LinearPlasticViscoelasticNormalSpecies::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")
106 {
107  return "LinearPlasticViscoelastic";
108 }

◆ getCohesionStiffness()

Mdouble LinearPlasticViscoelasticNormalSpecies::getCohesionStiffness ( ) const

Returns the cohesive stiffness of the linear plastic-viscoelastic normal force.

Returns
the cohesive stiffness of the linear plastic-viscoelastic normal force.
170 {
171  return cohesionStiffness_;
172 }

References cohesionStiffness_.

Referenced by LinearPlasticViscoelasticInteraction::computeNormalForce(), and mix().

◆ getCollisionTime()

Mdouble LinearPlasticViscoelasticNormalSpecies::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 If constant restitution is enabled, the collision time is mass-independent.

Todo:
should this use unloading stiffness?
318 {
319  if (getConstantRestitution()) mass = 1;
320 
321  logger.assert_debug(mass > 0, "getCollisionTime(%): mass should be positive",mass);
322 
323  Mdouble elasticContribution = loadingStiffness_ / (.5 * mass);
324  Mdouble elastoDissipativeContribution = elasticContribution - mathsFunc::square(dissipation_ / mass);
325 
326  logger.assert_debug(elastoDissipativeContribution > -1e-8 * elasticContribution,
327  "getCollisionTime(%): values for mass, stiffness and dissipation lead to an overdamped system: reduce dissipation",mass);
328 
329  return constants::pi / std::sqrt(elastoDissipativeContribution);
330 }

References dissipation_, BaseNormalForce::getConstantRestitution(), loadingStiffness_, logger, constants::pi, and mathsFunc::square().

Referenced by getRestitutionCoefficient().

◆ getDissipation()

Mdouble LinearPlasticViscoelasticNormalSpecies::getDissipation ( ) const

Allows the normal dissipation to be accessed.

Returns
the linear dissipation coefficient
249 {
250  return dissipation_;
251 }

References dissipation_.

Referenced by LinearPlasticViscoelasticInteraction::computeNormalForce(), mix(), and SlidingFrictionSpecies::setCollisionTimeAndNormalAndTangentialRestitutionCoefficient().

◆ getDoConstantUnloadingStiffness()

bool LinearPlasticViscoelasticNormalSpecies::getDoConstantUnloadingStiffness ( ) const
inline

◆ getLoadingStiffness()

Mdouble LinearPlasticViscoelasticNormalSpecies::getLoadingStiffness ( ) const

◆ getPenetrationDepthMax()

Mdouble LinearPlasticViscoelasticNormalSpecies::getPenetrationDepthMax ( ) const

Returns the maximum penetration depth of the linear plastic-viscoelastic normal force.

Returns
the maximum penetration depth of the linear plastic-viscoelastic normal force.
178 {
179  return penetrationDepthMax_;
180 }

References penetrationDepthMax_.

Referenced by LinearPlasticViscoelasticInteraction::computeNormalForce(), LinearPlasticViscoelasticInteraction::getUnloadingStiffness(), and mix().

◆ getRestitutionCoefficient()

Mdouble LinearPlasticViscoelasticNormalSpecies::getRestitutionCoefficient ( Mdouble  mass) const

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

333 {
334  if (getConstantRestitution()) mass = 1;
335 
336  return std::exp(-dissipation_ / mass * getCollisionTime(mass));
337 }
Mdouble getCollisionTime(Mdouble mass) const
Calculates collision time for two copies of a particle of given disp, k, mass.
Definition: LinearPlasticViscoelasticNormalSpecies.cc:317
Mdouble exp(Mdouble Exponent)
Definition: ExtendedMath.cc:84

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

◆ getUnloadingStiffnessMax()

Mdouble LinearPlasticViscoelasticNormalSpecies::getUnloadingStiffnessMax ( ) const

Returns the maximum unloading stiffness of the linear plastic-viscoelastic normal force.

Returns
the maximum unloading stiffness of the linear plastic-viscoelastic normal force.
162 {
163  return unloadingStiffnessMax_;
164 }

References unloadingStiffnessMax_.

Referenced by LinearPlasticViscoelasticInteraction::computeNormalForce(), LinearPlasticViscoelasticInteraction::getUnloadingStiffness(), and mix().

◆ mix()

void LinearPlasticViscoelasticNormalSpecies::mix ( LinearPlasticViscoelasticNormalSpecies S,
LinearPlasticViscoelasticNormalSpecies T 
)

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
117 {
123  logger.assert_debug(S->doConstantUnloadingStiffness_==T->doConstantUnloadingStiffness_,"you cannot mix species where doConstantUnloadingStiffness is not the same");
125 }
static Mdouble average(Mdouble a, Mdouble b)
Returns the harmonic mean of two variables.
Definition: BaseSpecies.cc:110
Mdouble getLoadingStiffness() const
Returns the loading stiffness of the linear plastic-viscoelastic normal force.
Definition: LinearPlasticViscoelasticNormalSpecies.cc:153
Mdouble getCohesionStiffness() const
Returns the cohesive stiffness of the linear plastic-viscoelastic normal force.
Definition: LinearPlasticViscoelasticNormalSpecies.cc:169
Mdouble getPenetrationDepthMax() const
Returns the maximum penetration depth of the linear plastic-viscoelastic normal force.
Definition: LinearPlasticViscoelasticNormalSpecies.cc:177
Mdouble getUnloadingStiffnessMax() const
Returns the maximum unloading stiffness of the linear plastic-viscoelastic normal force.
Definition: LinearPlasticViscoelasticNormalSpecies.cc:161
Mdouble getDissipation() const
Allows the normal dissipation to be accessed.
Definition: LinearPlasticViscoelasticNormalSpecies.cc:248

References BaseSpecies::average(), cohesionStiffness_, dissipation_, doConstantUnloadingStiffness_, getCohesionStiffness(), getDissipation(), getLoadingStiffness(), getPenetrationDepthMax(), getUnloadingStiffnessMax(), loadingStiffness_, logger, penetrationDepthMax_, and unloadingStiffnessMax_.

◆ read()

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

Reads the species properties from an input stream.

Parameters
[in]inputstream (typically the restart file)
92 {
93  std::string dummy;
94  is >> dummy >> loadingStiffness_;
95  is >> dummy >> unloadingStiffnessMax_;
96  helpers::readOptionalVariable<bool>(is,"doConstantUnloadingStiffness",doConstantUnloadingStiffness_);
97  is >> dummy >> cohesionStiffness_;
98  is >> dummy >> penetrationDepthMax_;
99  is >> dummy >> dissipation_;
100 }

References cohesionStiffness_, dissipation_, doConstantUnloadingStiffness_, loadingStiffness_, penetrationDepthMax_, and unloadingStiffnessMax_.

◆ setCohesionStiffness()

void LinearPlasticViscoelasticNormalSpecies::setCohesionStiffness ( Mdouble  cohesionStiffness)

Sets the cohesive stiffness of the linear plastic-viscoelastic normal force.

Parameters
[in]cohesionStiffnessthe cohesive stiffness of the linear plastic-viscoelastic normal force.
202 {
203  cohesionStiffness_ = cohesionStiffness;
204 }

References cohesionStiffness_.

Referenced by setPlasticParameters(), and ParticleParticleInteractionWithPlasticForces::setupInitialConditions().

◆ setCollisionTimeAndRestitutionCoefficient()

void LinearPlasticViscoelasticNormalSpecies::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)

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 mean of particle mass, \(\frac{2}{1/m1+1/m2}\) If constant restitution is enabled, the collision time and restitution is mass-independent.
Todo:
TW: check that the masses are described correctly here (m_eff or m_p?))
263 {
264  if (getConstantRestitution()) mass = 1;
265  if (eps == 0.0)
266  {
268  dissipation_ = std::sqrt(2.0 * mass * loadingStiffness_);
269  }
270  else
271  {
272  dissipation_ = -mass / tc * std::log(eps);
274  .5 * mass * (mathsFunc::square(constants::pi / tc) + mathsFunc::square(dissipation_ / mass));
275  }
277 }
Mdouble log(Mdouble Power)
Definition: ExtendedMath.cc:104

References dissipation_, BaseNormalForce::getConstantRestitution(), loadingStiffness_, mathsFunc::log(), constants::pi, mathsFunc::square(), and unloadingStiffnessMax_.

◆ setDissipation()

void LinearPlasticViscoelasticNormalSpecies::setDissipation ( Mdouble  dissipation)

Sets the linear dissipation coefficient of the linear plastic-viscoelastic normal force.

should be non-negative

Parameters
[in]thelinear dissipation coefficient of the linear plastic-viscoelastic normal force.
231 {
232  logger.assert_debug(dissipation >= 0, "setDissipation(%): dissipation should be non-negative",dissipation);
233  dissipation_ = dissipation;
234 }

References dissipation_, and logger.

Referenced by SlidingFrictionSpecies::setCollisionTimeAndNormalAndTangentialRestitutionCoefficient(), and setLoadingStiffnessAndDissipation().

◆ setDoConstantUnloadingStiffness()

void LinearPlasticViscoelasticNormalSpecies::setDoConstantUnloadingStiffness ( bool  doConstantUnloadingStiffness)
inline
171 {doConstantUnloadingStiffness_ = doConstantUnloadingStiffness;}

References doConstantUnloadingStiffness_.

◆ setLoadingStiffness()

void LinearPlasticViscoelasticNormalSpecies::setLoadingStiffness ( Mdouble  loadingStiffness)

Sets the loading stiffness of the linear plastic-viscoelastic normal force.

Parameters
[in]loadingStiffnessthe loading stiffness of the linear plastic-viscoelastic normal force.
186 {
187  loadingStiffness_ = loadingStiffness;
188 }

References loadingStiffness_.

Referenced by SlidingFrictionSpecies::setCollisionTimeAndNormalAndTangentialRestitutionCoefficient(), setLoadingStiffnessAndDissipation(), setPlasticParameters(), and ParticleParticleInteractionWithPlasticForces::setupInitialConditions().

◆ setLoadingStiffnessAndDissipation()

void LinearPlasticViscoelasticNormalSpecies::setLoadingStiffnessAndDissipation ( helpers::KAndDisp  new_)

Allows the spring and dissipation constants to be changed simultaneously.

Parameters
[in]ahelper struct containing both the loading stiffness and the disssipation coefficient.
240 {
241  setLoadingStiffness(new_.k);
242  setDissipation(new_.disp);
243 }
void setLoadingStiffness(Mdouble loadingStiffness)
Sets the loading stiffness of the linear plastic-viscoelastic normal force.
Definition: LinearPlasticViscoelasticNormalSpecies.cc:185
void setDissipation(Mdouble dissipation)
Sets the linear dissipation coefficient of the linear plastic-viscoelastic normal force.
Definition: LinearPlasticViscoelasticNormalSpecies.cc:230
Mdouble k
Definition: FormulaHelpers.h:39
Mdouble disp
Definition: FormulaHelpers.h:40

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

◆ setPenetrationDepthMax()

void LinearPlasticViscoelasticNormalSpecies::setPenetrationDepthMax ( Mdouble  penetrationDepthMax)

Sets the maximum penetration depth of the linear plastic-viscoelastic normal force.

Parameters
[in]penetrationDepthMaxthe maximum penetration depth of the linear plastic-viscoelastic normal force.
210 {
211  penetrationDepthMax_ = penetrationDepthMax;
212 }

References penetrationDepthMax_.

Referenced by setPlasticParameters(), and ParticleParticleInteractionWithPlasticForces::setupInitialConditions().

◆ setPlasticParameters()

void LinearPlasticViscoelasticNormalSpecies::setPlasticParameters ( Mdouble  loadingStiffness,
Mdouble  unloadingStiffnessMax,
Mdouble  cohesionStiffness,
Mdouble  penetrationDepthMax 
)

Sets all parameters of the linear plastic-viscoelastic normal force at once.

Parameters
[in]loadingStiffnessthe loading stiffness of the linear plastic-viscoelastic normal force.
[in]unloadingStiffnessMaxthe maximum unloading stiffness of the linear plastic-viscoelastic normal force.
[in]cohesionStiffnessthe cohesive stiffness of the linear plastic-viscoelastic normal force.
[in]penetrationDepthMaxthe maximum penetration depth of the linear plastic-viscoelastic normal force.
136 {
137  if (loadingStiffness <= 0 || unloadingStiffnessMax < loadingStiffness || cohesionStiffness < 0 ||
138  penetrationDepthMax < 0 || penetrationDepthMax > 1)
139  {
140  logger(ERROR,"Error: arguments of setPlasticParameters do not make sense\n %%%%%",
141  loadingStiffness <= 0, unloadingStiffnessMax < loadingStiffness, cohesionStiffness < 0,
142  penetrationDepthMax < 0, penetrationDepthMax > 1);
143  }
144  setLoadingStiffness(loadingStiffness);
145  setUnloadingStiffnessMax(unloadingStiffnessMax);
146  setCohesionStiffness(cohesionStiffness);
147  setPenetrationDepthMax(penetrationDepthMax);
148 }
@ ERROR
void setUnloadingStiffnessMax(Mdouble unloadingStiffnessMax)
Sets the maximum unloading stiffness of the linear plastic-viscoelastic normal force.
Definition: LinearPlasticViscoelasticNormalSpecies.cc:193
void setCohesionStiffness(Mdouble cohesionStiffness)
Sets the cohesive stiffness of the linear plastic-viscoelastic normal force.
Definition: LinearPlasticViscoelasticNormalSpecies.cc:201
void setPenetrationDepthMax(Mdouble penetrationDepthMax)
Sets the maximum penetration depth of the linear plastic-viscoelastic normal force.
Definition: LinearPlasticViscoelasticNormalSpecies.cc:209

References ERROR, logger, setCohesionStiffness(), setLoadingStiffness(), setPenetrationDepthMax(), and setUnloadingStiffnessMax().

◆ setRestitutionCoefficient()

void LinearPlasticViscoelasticNormalSpecies::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}\) If constant restitution is enabled, the restitution coefficient is mass-independent.
301 {
302  if (getConstantRestitution()) mass = 1;
303  if (eps == 0.0) {
304  dissipation_ = std::sqrt(2.0 * mass * getLoadingStiffness());
305  } else {
306  const Mdouble logEps = log(eps);
307  dissipation_ = -std::sqrt(2.0 * mass * getLoadingStiffness()
308  / (constants::sqr_pi + mathsFunc::square(logEps))) * logEps;
309  }
310 }
const Mdouble sqr_pi
Definition: ExtendedMath.h:47

References dissipation_, BaseNormalForce::getConstantRestitution(), getLoadingStiffness(), mathsFunc::log(), constants::sqr_pi, and mathsFunc::square().

Referenced by setStiffnessAndRestitutionCoefficient().

◆ setStiffnessAndRestitutionCoefficient()

void LinearPlasticViscoelasticNormalSpecies::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]masseffective particle mass, \(\frac{2}{1/m1+1/m2}\) If constant restitution is enabled, the restitution coefficient is mass-independent.
287 {
288  if (getConstantRestitution()) mass = 1;
289  loadingStiffness_ = stiffness;
290  setRestitutionCoefficient(eps, mass);
291 }
void setRestitutionCoefficient(double eps, Mdouble mass)
Sets disp to obtain a restitution coefficient eps for a collision of two particles of mass m.
Definition: LinearPlasticViscoelasticNormalSpecies.cc:300

References BaseNormalForce::getConstantRestitution(), loadingStiffness_, and setRestitutionCoefficient().

◆ setUnloadingStiffnessMax()

void LinearPlasticViscoelasticNormalSpecies::setUnloadingStiffnessMax ( Mdouble  unloadingStiffnessMax)

Sets the maximum unloading stiffness of the linear plastic-viscoelastic normal force.

Parameters
[in]unloadingStiffnessMaxthe maximum unloading stiffness of the linear plastic-viscoelastic normal force.
194 {
195  unloadingStiffnessMax_ = unloadingStiffnessMax;
196 }

References unloadingStiffnessMax_.

Referenced by setPlasticParameters(), and ParticleParticleInteractionWithPlasticForces::setupInitialConditions().

◆ write()

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

Writes the species properties to an output stream.

Parameters
[out]outputstream (typically the restart file)
79 {
80  os << " loadingStiffness " << loadingStiffness_;
81  os << " maxUnloadingStiffness " << unloadingStiffnessMax_;
82  if (doConstantUnloadingStiffness_) os << " doConstantUnloadingStiffness " << doConstantUnloadingStiffness_;
83  os << " cohesionStiffness " << cohesionStiffness_;
84  os << " maxPenetration " << penetrationDepthMax_;
85  os << " dissipation " << dissipation_;
86 }

References cohesionStiffness_, dissipation_, doConstantUnloadingStiffness_, loadingStiffness_, penetrationDepthMax_, and unloadingStiffnessMax_.

Member Data Documentation

◆ cohesionStiffness_

Mdouble LinearPlasticViscoelasticNormalSpecies::cohesionStiffness_
private

the adhesive spring constant (k^c) for plastic deformations

Referenced by computeBondNumberMax(), getCohesionStiffness(), LinearPlasticViscoelasticNormalSpecies(), mix(), read(), setCohesionStiffness(), and write().

◆ dissipation_

◆ doConstantUnloadingStiffness_

bool LinearPlasticViscoelasticNormalSpecies::doConstantUnloadingStiffness_ = false
private

◆ loadingStiffness_

◆ penetrationDepthMax_

Mdouble LinearPlasticViscoelasticNormalSpecies::penetrationDepthMax_
private

the depth (relative to the normalized radius) at which k_2^max is used (phi_f)

Referenced by computeBondNumberMax(), getPenetrationDepthMax(), LinearPlasticViscoelasticNormalSpecies(), mix(), read(), setPenetrationDepthMax(), and write().

◆ unloadingStiffnessMax_

Mdouble LinearPlasticViscoelasticNormalSpecies::unloadingStiffnessMax_
private

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