SinterLinNormalSpecies Class Reference

SinterLinNormalSpecies contains the parameters used to describe a plastic-cohesive normal force (Stefan Ludings plastic-cohesive force model) based on three different sintering mechanisms. More...

#include <SinterLinNormalSpecies.h>

+ Inheritance diagram for SinterLinNormalSpecies:

Public Types

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

Public Member Functions

 SinterLinNormalSpecies ()
 The default constructor. More...
 
 SinterLinNormalSpecies (const SinterLinNormalSpecies &p)
 The default copy constructor. More...
 
 ~SinterLinNormalSpecies ()
 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 (SinterLinNormalSpecies *S, SinterLinNormalSpecies *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...
 
void setSinterAdhesion (Mdouble sinterAdhesion)
 Sets sinterAdhesion_. More...
 
Mdouble getSinterAdhesion () const
 Accesses sinterAdhesion_. More...
 
void setSinterRate (Mdouble sinterRate)
 Sets sinterRate_. More...
 
Mdouble getSinterRate () const
 Accesses sinterRate_. More...
 
void setComplianceZero (Mdouble complianceZero)
 sets the instantaneous compliance (compliance zero) More...
 
Mdouble getComplianceZero () const
 Accesses the instantaneous compliance (compliance zero) More...
 
void setSurfTension (Mdouble complianceZero)
 sets the surface tension. More...
 
Mdouble getSurfTension () const
 accesses the surface tension. More...
 
void setFluidity (Mdouble complianceOne)
 sets the fluidity (inverse of viscosity). More...
 
Mdouble getFluidity () const
 accesses the fluidity (inverse of viscosity). More...
 
void setSeparationDis (Mdouble separationDis)
 sets the critical separation distance. This mimics inter-surface forces. More...
 
Mdouble getSeparationDis () const
 accesses the critical separation distance. This mimics inter-surface forces. More...
 
void setSinterType (SINTER_APPROACH sinterType)
 sets the sintering approach selected. More...
 
SINTER_APPROACH getSinterType () const
 
- 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...
 
Mdouble sinterAdhesion_
 
Mdouble sinterRate_
 
Mdouble complianceZero_
 
Mdouble surfTension_
 
Mdouble fluidity
 
Mdouble separationDis_
 
SINTER_APPROACH sinterType_
 Type of sintering. More...
 

Detailed Description

SinterLinNormalSpecies contains the parameters used to describe a plastic-cohesive normal force (Stefan Ludings plastic-cohesive force model) based on three different sintering mechanisms.

Ref: Visco-elastic sintering kinetics in virgin and aged polymer powders. Powder Technology, 2020.

Member Typedef Documentation

◆ InteractionType

The correct Interaction type for this FrictionForceSpecies.

Constructor & Destructor Documentation

◆ SinterLinNormalSpecies() [1/2]

SinterLinNormalSpecies::SinterLinNormalSpecies ( )

The default constructor.

37  : BaseNormalForce()
38 {
39  loadingStiffness_ = 0.0;
41  cohesionStiffness_ = 0.0;
43  dissipation_ = 0.0;
44  sinterAdhesion_ = 0.0;
45  sinterRate_ = 0.0;
46  complianceZero_ = 0.0;
47  surfTension_ = 0.0;
48  fluidity = 0.0;
49  separationDis_ = 0.0;
51 
52 #ifdef DEBUG_CONSTRUCTOR
53  std::cout<<"SinterLinNormalSpecies::SinterLinNormalSpecies() finished"<<std::endl;
54 #endif
55 }
BaseNormalForce()
Definition: BaseNormalForce.h:35
Mdouble separationDis_
Definition: SinterLinNormalSpecies.h:290
SINTER_APPROACH sinterType_
Type of sintering.
Definition: SinterLinNormalSpecies.h:293
Mdouble cohesionStiffness_
the adhesive spring constant (k^c) for plastic deformations
Definition: SinterLinNormalSpecies.h:253
Mdouble complianceZero_
Definition: SinterLinNormalSpecies.h:275
Mdouble fluidity
Definition: SinterLinNormalSpecies.h:285
Mdouble sinterRate_
Definition: SinterLinNormalSpecies.h:269
Mdouble penetrationDepthMax_
the depth (relative to the normalized radius) at which k_2^max is used (phi_f)
Definition: SinterLinNormalSpecies.h:256
Mdouble unloadingStiffnessMax_
the maximum elastic constant (k_2^max) for plastic deformations
Definition: SinterLinNormalSpecies.h:250
Mdouble sinterAdhesion_
Definition: SinterLinNormalSpecies.h:264
Mdouble surfTension_
Definition: SinterLinNormalSpecies.h:280
Mdouble dissipation_
linear dissipation coefficient
Definition: SinterLinNormalSpecies.h:259
Mdouble loadingStiffness_
(normal) spring constant (k_1)
Definition: SinterLinNormalSpecies.h:247

References cohesionStiffness_, complianceZero_, dissipation_, fluidity, FRENKEL, loadingStiffness_, penetrationDepthMax_, separationDis_, sinterAdhesion_, sinterRate_, sinterType_, surfTension_, and unloadingStiffnessMax_.

◆ SinterLinNormalSpecies() [2/2]

SinterLinNormalSpecies::SinterLinNormalSpecies ( const SinterLinNormalSpecies p)

The default copy constructor.

Parameters
[in]thespecies that is copied
61  : BaseNormalForce(p)
62 {
72  fluidity = p.fluidity;
75 
76 #ifdef DEBUG_CONSTRUCTOR
77  std::cout<<"SinterLinNormalSpecies::SinterLinNormalSpecies(const SinterLinNormalSpecies &p) finished"<<std::endl;
78 #endif
79 }

References cohesionStiffness_, complianceZero_, dissipation_, fluidity, loadingStiffness_, penetrationDepthMax_, separationDis_, sinterAdhesion_, sinterRate_, sinterType_, surfTension_, and unloadingStiffnessMax_.

◆ ~SinterLinNormalSpecies()

SinterLinNormalSpecies::~SinterLinNormalSpecies ( )

The default destructor.

82 {
83 #ifdef DEBUG_DESTRUCTOR
84  std::cout<<"SinterNormalSpecies::~SinterNormalSpecies() finished"<<std::endl;
85 #endif
86 }

Member Function Documentation

◆ computeTimeStep()

Mdouble SinterLinNormalSpecies::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]massthe optimal time step is computed to resolve a collision of two particles of this mass.
250 {
251  if (getConstantRestitution()) mass = 1;
252  return 0.02 * constants::pi /
253  std::sqrt(unloadingStiffnessMax_ / (.5 * mass) - mathsFunc::square(dissipation_ / mass));
254 }
bool getConstantRestitution() const
Definition: BaseNormalForce.h:46
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 SinterLinNormalSpecies::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")
134 {
135  return "SinterLin";
136 }

◆ getCohesionStiffness()

Mdouble SinterLinNormalSpecies::getCohesionStiffness ( ) const

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

Returns
the cohesive stiffness of the linear plastic-viscoelastic normal force.
200 {
201  return cohesionStiffness_;
202 }

References cohesionStiffness_.

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

◆ getCollisionTime()

Mdouble SinterLinNormalSpecies::getCollisionTime ( Mdouble  mass) const

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

457 {
458  if (mass <= 0)
459  {
460  std::cerr << "Error in getCollisionTime(" << mass
461  << ") mass is not set or has an unexpected value, (getCollisionTime(" << mass << "))" << std::endl;
462  exit(-1);
463  }
464  if (loadingStiffness_ <= 0)
465  {
466  std::cerr << "Error in getCollisionTime(" << mass << ") stiffness=" << loadingStiffness_
467  << " is not set or has an unexpected value, (getCollisionTime(" << mass << "), with stiffness="
468  << loadingStiffness_ << ")" << std::endl;
469  exit(-1);
470  }
471  if (dissipation_ < 0)
472  {
473  std::cerr << "Error in getCollisionTime(" << mass << ") dissipation=" << dissipation_
474  << " is not set or has an unexpected value, (getCollisionTime(" << mass << "), with dissipation="
475  << dissipation_ << ")" << std::endl;
476  exit(-1);
477  }
478  Mdouble tosqrt = loadingStiffness_ / (.5 * mass) - mathsFunc::square(dissipation_ / mass);
479  if (tosqrt <= -1e-8 * loadingStiffness_ / (.5 * mass))
480  {
481  std::cerr << "Error in getCollisionTime(" << mass
482  << ") values for mass, stiffness and dissipation would lead to an overdamped system, (getCollisionTime("
483  << mass << "), with stiffness=" << loadingStiffness_ << " and dissipation=" << dissipation_ << ")"
484  << std::endl;
485  exit(-1);
486  }
487  return constants::pi / std::sqrt(tosqrt);
488 }
double Mdouble
Definition: GeneralDefine.h:34

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

Referenced by getRestitutionCoefficient().

◆ getComplianceZero()

Mdouble SinterLinNormalSpecies::getComplianceZero ( ) const

Accesses the instantaneous compliance (compliance zero)

372 {
373  return complianceZero_;
374 }

References complianceZero_.

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

◆ getDissipation()

Mdouble SinterLinNormalSpecies::getDissipation ( ) const

Allows the normal dissipation to be accessed.

Returns
the linear dissipation coefficient
351 {
352  return dissipation_;
353 }

References dissipation_.

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

◆ getFluidity()

Mdouble SinterLinNormalSpecies::getFluidity ( ) const

accesses the fluidity (inverse of viscosity).

382 {
383  return fluidity;
384 }

References fluidity.

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

◆ getLoadingStiffness()

Mdouble SinterLinNormalSpecies::getLoadingStiffness ( ) const

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

Returns
the loading stiffness of the linear plastic-viscoelastic normal force.
184 {
185  return loadingStiffness_;
186 }

References loadingStiffness_.

Referenced by SinterLinInteraction::computeNormalForce(), SinterLinInteraction::getElasticEnergy(), SinterLinInteraction::getUnloadingStiffness(), mix(), and setRestitutionCoefficient().

◆ getPenetrationDepthMax()

Mdouble SinterLinNormalSpecies::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.
208 {
209  return penetrationDepthMax_;
210 }

References penetrationDepthMax_.

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

◆ getRestitutionCoefficient()

Mdouble SinterLinNormalSpecies::getRestitutionCoefficient ( Mdouble  mass) const

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

491 {
492  if (getConstantRestitution()) mass = 1;
493 
494  return std::exp(-dissipation_ / mass * getCollisionTime(mass));
495 }
Mdouble getCollisionTime(Mdouble mass) const
Calculates collision time for two copies of a particle of given disp, k, mass.
Definition: SinterLinNormalSpecies.cc:456
Mdouble exp(Mdouble Exponent)
Definition: ExtendedMath.cc:84

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

◆ getSeparationDis()

Mdouble SinterLinNormalSpecies::getSeparationDis ( ) const

accesses the critical separation distance. This mimics inter-surface forces.

387 {
388  return separationDis_;
389 }

References separationDis_.

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

◆ getSinterAdhesion()

Mdouble SinterLinNormalSpecies::getSinterAdhesion ( ) const

Accesses sinterAdhesion_.

Returns
value of sinterAdhesion_
359 {
360  return sinterAdhesion_;
361 }

References sinterAdhesion_.

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

◆ getSinterRate()

Mdouble SinterLinNormalSpecies::getSinterRate ( ) const

Accesses sinterRate_.

Returns
value of sinterAdhesion_
367 {
368  return sinterRate_;
369 }

References sinterRate_.

Referenced by SinterLinInteraction::computeNormalForce().

◆ getSinterType()

SINTER_APPROACH SinterLinNormalSpecies::getSinterType ( ) const
393 {
394  return sinterType_;
395 }

References sinterType_.

Referenced by SinterLinInteraction::computeNormalForce().

◆ getSurfTension()

Mdouble SinterLinNormalSpecies::getSurfTension ( ) const

accesses the surface tension.

377 {
378  return surfTension_;
379 }

References surfTension_.

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

◆ getUnloadingStiffnessMax()

Mdouble SinterLinNormalSpecies::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.
192 {
193  return unloadingStiffnessMax_;
194 }

References unloadingStiffnessMax_.

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

◆ mix()

void SinterLinNormalSpecies::mix ( SinterLinNormalSpecies S,
SinterLinNormalSpecies 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
145 {
156 
157 }
static Mdouble average(Mdouble a, Mdouble b)
Returns the harmonic mean of two variables.
Definition: BaseSpecies.cc:110
Mdouble getUnloadingStiffnessMax() const
Returns the maximum unloading stiffness of the linear plastic-viscoelastic normal force.
Definition: SinterLinNormalSpecies.cc:191
Mdouble getSeparationDis() const
accesses the critical separation distance. This mimics inter-surface forces.
Definition: SinterLinNormalSpecies.cc:386
Mdouble getCohesionStiffness() const
Returns the cohesive stiffness of the linear plastic-viscoelastic normal force.
Definition: SinterLinNormalSpecies.cc:199
Mdouble getPenetrationDepthMax() const
Returns the maximum penetration depth of the linear plastic-viscoelastic normal force.
Definition: SinterLinNormalSpecies.cc:207
Mdouble getSinterAdhesion() const
Accesses sinterAdhesion_.
Definition: SinterLinNormalSpecies.cc:358
Mdouble getLoadingStiffness() const
Returns the loading stiffness of the linear plastic-viscoelastic normal force.
Definition: SinterLinNormalSpecies.cc:183
Mdouble getComplianceZero() const
Accesses the instantaneous compliance (compliance zero)
Definition: SinterLinNormalSpecies.cc:371
Mdouble getSurfTension() const
accesses the surface tension.
Definition: SinterLinNormalSpecies.cc:376
Mdouble getFluidity() const
accesses the fluidity (inverse of viscosity).
Definition: SinterLinNormalSpecies.cc:381
Mdouble getDissipation() const
Allows the normal dissipation to be accessed.
Definition: SinterLinNormalSpecies.cc:350

References BaseSpecies::average(), cohesionStiffness_, complianceZero_, dissipation_, fluidity, getCohesionStiffness(), getComplianceZero(), getDissipation(), getFluidity(), getLoadingStiffness(), getPenetrationDepthMax(), getSeparationDis(), getSinterAdhesion(), getSurfTension(), getUnloadingStiffnessMax(), loadingStiffness_, penetrationDepthMax_, separationDis_, sinterAdhesion_, surfTension_, and unloadingStiffnessMax_.

◆ read()

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

Reads the species properties from an input stream.

Parameters
[in]isinput stream (typically the restart file)
111 {
112  std::string dummy;
113  is >> dummy >> loadingStiffness_;
114  is >> dummy >> unloadingStiffnessMax_;
115  is >> dummy >> cohesionStiffness_;
116  is >> dummy >> penetrationDepthMax_;
117  is >> dummy >> dissipation_;
118  is >> dummy >> sinterAdhesion_;
119  is >> dummy >> complianceZero_;
120  is >> dummy >> surfTension_;
121  is >> dummy >> sinterRate_;
122  is >> dummy >> fluidity;
123  is >> dummy >> separationDis_;
124  unsigned type;
125  is >> dummy >> type;
126 
127  sinterType_ = (SINTER_APPROACH) type;
128 }
SINTER_APPROACH
Definition: SinterLinNormalSpecies.h:40

References cohesionStiffness_, complianceZero_, dissipation_, fluidity, loadingStiffness_, penetrationDepthMax_, separationDis_, sinterAdhesion_, sinterRate_, sinterType_, surfTension_, and unloadingStiffnessMax_.

◆ setCohesionStiffness()

void SinterLinNormalSpecies::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.
232 {
233  cohesionStiffness_ = cohesionStiffness;
234 }

References cohesionStiffness_.

Referenced by setPlasticParameters().

◆ setCollisionTimeAndRestitutionCoefficient()

void SinterLinNormalSpecies::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]masseffective particle mass, \(\frac{2}{1/m1+1/m2}\)
Todo:
TW: check that the masses are described correctly here (m_eff or m_p?))
407 {
408  if (eps == 0.0)
409  {
411  dissipation_ = std::sqrt(2.0 * mass * loadingStiffness_);
412  }
413  else
414  {
415  dissipation_ = -mass / tc * std::log(eps);
417  .5 * mass * (mathsFunc::square(constants::pi / tc) + mathsFunc::square(dissipation_ / mass));
418  }
420 }
Mdouble log(Mdouble Power)
Definition: ExtendedMath.cc:104

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

◆ setComplianceZero()

void SinterLinNormalSpecies::setComplianceZero ( Mdouble  complianceZero)

sets the instantaneous compliance (compliance zero)

297 {
298  if (complianceZero < 0)
299  {
300  logger(ERROR, "complianceZero(%)", complianceZero);
301  }
302  complianceZero_ = complianceZero;
303 }
Logger< MERCURYDPM_LOGLEVEL > logger("MercuryKernel")
Definition of different loggers with certain modules. A user can define its own custom logger here.
@ ERROR

References complianceZero_, ERROR, and logger.

◆ setDissipation()

void SinterLinNormalSpecies::setDissipation ( Mdouble  dissipation)

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

should be non-negative

Parameters
[in]dissipationthe linear dissipation coefficient of the linear plastic-viscoelastic normal force.
261 {
262  if (dissipation < 0)
263  {
264  logger(ERROR, "setDissipation(%)", dissipation);
265  exit(-1);
266  }
267  dissipation_ = dissipation;
268 }

References dissipation_, ERROR, and logger.

Referenced by setLoadingStiffnessAndDissipation().

◆ setFluidity()

void SinterLinNormalSpecies::setFluidity ( Mdouble  complianceOne)

sets the fluidity (inverse of viscosity).

315 {
316  if (complianceOne < 0)
317  {
318  logger(ERROR, "ComplianceOne(%)", complianceOne);
319  }
320  fluidity = complianceOne;
321 }

References ERROR, fluidity, and logger.

◆ setLoadingStiffness()

void SinterLinNormalSpecies::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.
216 {
217  loadingStiffness_ = loadingStiffness;
218 }

References loadingStiffness_.

Referenced by setLoadingStiffnessAndDissipation(), and setPlasticParameters().

◆ setLoadingStiffnessAndDissipation()

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

Allows the spring and dissipation constants to be changed simultaneously.

272 {
273  setLoadingStiffness(new_.k);
274  setDissipation(new_.disp);
275 }
void setDissipation(Mdouble dissipation)
Sets the linear dissipation coefficient of the linear plastic-viscoelastic normal force.
Definition: SinterLinNormalSpecies.cc:260
void setLoadingStiffness(Mdouble loadingStiffness)
Sets the loading stiffness of the linear plastic-viscoelastic normal force.
Definition: SinterLinNormalSpecies.cc:215
Mdouble k
Definition: FormulaHelpers.h:39
Mdouble disp
Definition: FormulaHelpers.h:40

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

◆ setPenetrationDepthMax()

void SinterLinNormalSpecies::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.
240 {
241  penetrationDepthMax_ = penetrationDepthMax;
242 }

References penetrationDepthMax_.

Referenced by setPlasticParameters().

◆ setPlasticParameters()

void SinterLinNormalSpecies::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.
167 {
168  if (loadingStiffness <= 0 || unloadingStiffnessMax < loadingStiffness || cohesionStiffness < 0 ||
169  penetrationDepthMax < 0)
170  {
171  std::cerr << "Error: arguments of setPlasticParameters do not make sense" << std::endl;
172  exit(-1);
173  }
174  setLoadingStiffness(loadingStiffness);
175  setUnloadingStiffnessMax(unloadingStiffnessMax);
176  setCohesionStiffness(cohesionStiffness);
177  setPenetrationDepthMax(penetrationDepthMax);
178 }
void setCohesionStiffness(Mdouble cohesionStiffness)
Sets the cohesive stiffness of the linear plastic-viscoelastic normal force.
Definition: SinterLinNormalSpecies.cc:231
void setUnloadingStiffnessMax(Mdouble unloadingStiffnessMax)
Sets the maximum unloading stiffness of the linear plastic-viscoelastic normal force.
Definition: SinterLinNormalSpecies.cc:223
void setPenetrationDepthMax(Mdouble penetrationDepthMax)
Sets the maximum penetration depth of the linear plastic-viscoelastic normal force.
Definition: SinterLinNormalSpecies.cc:239

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

◆ setRestitutionCoefficient()

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

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

444 {
445  if (getConstantRestitution()) mass = 1;
446  if (eps == 0.0) {
447  dissipation_ = std::sqrt(2.0 * mass * getLoadingStiffness());
448  } else {
449  const Mdouble logEps = log(eps);
450  dissipation_ = -std::sqrt(2.0 * mass * getLoadingStiffness()
451  / (constants::sqr_pi + mathsFunc::square(logEps))) * logEps;
452  }
453 }
const Mdouble sqr_pi
Definition: ExtendedMath.h:47

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

◆ setSeparationDis()

void SinterLinNormalSpecies::setSeparationDis ( Mdouble  separationDis)

sets the critical separation distance. This mimics inter-surface forces.

324 {
325  if (separationDis < 0)
326  {
327  logger(ERROR, "SeparationDis(%)", separationDis);
328  }
329  separationDis_ = separationDis;
330 }

References ERROR, logger, and separationDis_.

◆ setSinterAdhesion()

void SinterLinNormalSpecies::setSinterAdhesion ( Mdouble  sinterAdhesion)

Sets sinterAdhesion_.

278 {
279  if (sinterAdhesion < 0)
280  {
281  logger(ERROR, "setSinterAdhesion(%)", sinterAdhesion);
282  exit(-1);
283  }
284  sinterAdhesion_ = sinterAdhesion;
285 }

References ERROR, logger, and sinterAdhesion_.

◆ setSinterRate()

void SinterLinNormalSpecies::setSinterRate ( Mdouble  sinterRate)

Sets sinterRate_.

288 {
289  if (sinterRate < 0)
290  {
291  logger(ERROR, "setSinterRate(%)", sinterRate);
292  }
293  sinterRate_ = sinterRate;
294 }

References ERROR, logger, and sinterRate_.

◆ setSinterType()

void SinterLinNormalSpecies::setSinterType ( SINTER_APPROACH  sinterType)

sets the sintering approach selected.

333 {
334  sinterType_ = sinterType;
336  logger(INFO, "Sintertype set to Frenkel");
338  logger(INFO, "Sintertype set to Three different mechanisms");
339  else
340  logger(ERROR, "Sintertype not understood");
341 }
@ INFO

References ERROR, FRENKEL, INFO, logger, sinterType_, and VISCOELASTIC_CONTACT.

◆ setStiffnessAndRestitutionCoefficient()

void SinterLinNormalSpecies::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}\)
429 {
430  loadingStiffness_ = stiffness;
431  if (eps == 0.0)
432  {
433  dissipation_ = std::sqrt(2.0 * mass * stiffness);
434  }
435  else
436  {
437  dissipation_ =
438  -std::sqrt(2.0 * mass * stiffness / (constants::sqr_pi + mathsFunc::square(log(eps)))) * log(eps);
439  }
441 }

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

◆ setSurfTension()

void SinterLinNormalSpecies::setSurfTension ( Mdouble  complianceZero)

sets the surface tension.

306 {
307  if (surfTension < 0)
308  {
309  logger(ERROR, "SurfTension(%)", surfTension);
310  }
311  surfTension_ = surfTension;
312 }

References ERROR, logger, and surfTension_.

◆ setUnloadingStiffnessMax()

void SinterLinNormalSpecies::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.
224 {
225  unloadingStiffnessMax_ = unloadingStiffnessMax;
226 }

References unloadingStiffnessMax_.

Referenced by setPlasticParameters().

◆ write()

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

Writes the species properties to an output stream.

Parameters
[out]outputstream (typically the restart file)
92 {
93  os << " loadingStiffness " << loadingStiffness_;
94  os << " maxUnloadingStiffness " << unloadingStiffnessMax_;
95  os << " cohesionStiffness " << cohesionStiffness_;
96  os << " maxPenetration " << penetrationDepthMax_;
97  os << " dissipation " << dissipation_;
98  os << " sinterAdhesion " << sinterAdhesion_;
99  os << " sinterRate " << sinterRate_;
100  os << " complianceZero " << complianceZero_;
101  os << " surfTension " << surfTension_;
102  os << " complianceOne " << fluidity;
103  os << " separationDis " << separationDis_;
104  os << " sinterType " << (unsigned) sinterType_;
105 }

References cohesionStiffness_, complianceZero_, dissipation_, fluidity, loadingStiffness_, penetrationDepthMax_, separationDis_, sinterAdhesion_, sinterRate_, sinterType_, surfTension_, and unloadingStiffnessMax_.

Member Data Documentation

◆ cohesionStiffness_

Mdouble SinterLinNormalSpecies::cohesionStiffness_
private

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

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

◆ complianceZero_

Mdouble SinterLinNormalSpecies::complianceZero_
private

Material property. C0 = (1-\nu^2)/E, E: young's Modulus, nu poisson ratio.

Referenced by getComplianceZero(), mix(), read(), setComplianceZero(), SinterLinNormalSpecies(), and write().

◆ dissipation_

◆ fluidity

Mdouble SinterLinNormalSpecies::fluidity
private

Material property, which represents the inverse of viscosity

Referenced by getFluidity(), mix(), read(), setFluidity(), SinterLinNormalSpecies(), and write().

◆ loadingStiffness_

◆ penetrationDepthMax_

Mdouble SinterLinNormalSpecies::penetrationDepthMax_
private

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

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

◆ separationDis_

Mdouble SinterLinNormalSpecies::separationDis_
private

Critical separation distance, to describe inter-surface forces.

Referenced by getSeparationDis(), mix(), read(), setSeparationDis(), SinterLinNormalSpecies(), and write().

◆ sinterAdhesion_

Mdouble SinterLinNormalSpecies::sinterAdhesion_
private

Determines sinter adhesion force fa=sinterAdhesion_*radius in sinter rate: d(delta0)/dt = (fa+fep)/nu

Referenced by getSinterAdhesion(), mix(), read(), setSinterAdhesion(), SinterLinNormalSpecies(), and write().

◆ sinterRate_

Mdouble SinterLinNormalSpecies::sinterRate_
private

Determines sinter rate: d(delta0)/dt = (fa+fep)/nu

Referenced by getSinterRate(), read(), setSinterRate(), SinterLinNormalSpecies(), and write().

◆ sinterType_

SINTER_APPROACH SinterLinNormalSpecies::sinterType_
private

Type of sintering.

Referenced by getSinterType(), read(), setSinterType(), SinterLinNormalSpecies(), and write().

◆ surfTension_

Mdouble SinterLinNormalSpecies::surfTension_
private

Material property.

Referenced by getSurfTension(), mix(), read(), setSurfTension(), SinterLinNormalSpecies(), and write().

◆ unloadingStiffnessMax_

Mdouble SinterLinNormalSpecies::unloadingStiffnessMax_
private

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