MercuryDPM  Alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
SinterNormalSpecies Class Reference

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

#include <SinterNormalSpecies.h>

+ Inheritance diagram for SinterNormalSpecies:

Public Types

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

Public Member Functions

 SinterNormalSpecies ()
 The default constructor. More...
 
 SinterNormalSpecies (const SinterNormalSpecies &p)
 The default copy constructor. More...
 
virtual ~SinterNormalSpecies ()
 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 (SinterNormalSpecies *const S, SinterNormalSpecies *const T)
 creates default values for mixed species More...
 
void 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. 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...
 
Mdouble getCollisionTime (Mdouble mass)
 Calculates collision time for two copies of a particle 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...
 
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 setInverseSinterViscosity (Mdouble inverseSinterViscosity)
 Sets inverseSinterViscosity_. More...
 
Mdouble getInverseSinterViscosity () const
 Accesses inverseSinterViscosity_. More...
 
MERCURY_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...
 
void setSinterForceAndTime (Mdouble adhesionForce, Mdouble sinterTime, Mdouble radius)
 
void setParhamiMcKeeping (Mdouble alpha, Mdouble beta, Mdouble atomicVolume, Mdouble surfaceEnergy, Mdouble thicknessDiffusion, Mdouble activationEnergy, Mdouble temperature)
 Sets the sinterAdhesion_ and inverseSinterViscosity_ based on the Parhami-McKeeping parameters. More...
 
void setSinterRate (Mdouble sinterRate)
 Sets sinterRate_. More...
 
void setSinterType (SINTERTYPE sinterType)
 Sets sinterRate_. More...
 
Mdouble getSinterRate () const
 Accesses sinterRate_. More...
 
SINTERTYPE getSinterType () const
 
std::function< double(double
temperature)> 
getTemperatureDependentSinterRate () const
 
double getTemperatureDependentSinterRate (double temperature) const
 
void setTemperatureDependentSinterRate (std::function< double(double temperature)> temperatureDependentSinterRate)
 
- Public Member Functions inherited from BaseSpecies
 BaseSpecies ()
 The default constructor. More...
 
 BaseSpecies (const BaseSpecies &p)
 The copy constructor. More...
 
virtual ~BaseSpecies ()
 The default destructor. More...
 
virtual BaseSpeciescopy () const =0
 Creates a deep copy of the object from which it is called. More...
 
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...
 
Mdouble average (Mdouble a, Mdouble b)
 defines the average of two variables by the harmonic mean. More...
 
virtual void mixAll (BaseSpecies *const S, BaseSpecies *const T)=0
 creates default values for mixed species More...
 
virtual Mdouble getInteractionDistance () const =0
 returns the largest separation distance at which adhesive short-range forces can occur. 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, Mdouble timeStamp) const =0
 returns new Interaction object. More...
 
- Public Member Functions inherited from BaseObject
 BaseObject ()
 Default constructor. More...
 
 BaseObject (const BaseObject &p)
 Copy constructor, copies all the objects BaseObject contains. More...
 
virtual ~BaseObject ()
 virtual destructor More...
 
virtual std::string getName () const =0
 A purely virtual function. More...
 
virtual void moveInHandler (const unsigned int index)
 Except that it is virtual, it does the same thing as setIndex() does. More...
 
void setIndex (const unsigned int index)
 Allows one to assign an index to an object in the handler/container. More...
 
void setId (const unsigned int 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...
 

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 inverseSinterViscosity_
 
Mdouble sinterRate_
 
SINTERTYPE sinterType_
 sinterType options determin how the rate of sintering, d(equilibriumOverlap)/dt is computed PARHAMI_MCKEEPING: sinter rate given by (sinterAdhesion+normalForce)/sinterViscosity CONSTANT_RATE: sinter rate given by sinterRate More...
 
std::function< double(double
temperature)> 
temperatureDependentSinterRate_ = [this] (double temperature) {return sinterRate_;}
 

Detailed Description

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

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

Definition at line 42 of file SinterNormalSpecies.h.

Member Typedef Documentation

Constructor & Destructor Documentation

SinterNormalSpecies::SinterNormalSpecies ( )

The default constructor.

Definition at line 33 of file SinterNormalSpecies.cc.

References cohesionStiffness_, dissipation_, inverseSinterViscosity_, loadingStiffness_, PARHAMI_MCKEEPING, penetrationDepthMax_, sinterAdhesion_, sinterRate_, sinterType_, and unloadingStiffnessMax_.

34 {
36  loadingStiffness_ = 0.0;
38  cohesionStiffness_ = 0.0;
40  dissipation_ = 0.0;
41  sinterAdhesion_ = 0.0;
42  sinterRate_ = 0.0;
44 #ifdef DEBUG_CONSTRUCTOR
45  std::cout<<"SinterNormalSpecies::SinterNormalSpecies() finished"<<std::endl;
46 #endif
47 }
Mdouble unloadingStiffnessMax_
the maximum elastic constant (k_2^max) for plastic deformations
SINTERTYPE sinterType_
sinterType options determin how the rate of sintering, d(equilibriumOverlap)/dt is computed PARHAMI_M...
Mdouble loadingStiffness_
(normal) spring constant (k_1)
Mdouble cohesionStiffness_
the adhesive spring constant (k^c) for plastic deformations
Mdouble penetrationDepthMax_
the depth (relative to the normalized radius) at which k_2^max is used (phi_f)
Mdouble dissipation_
linear dissipation coefficient
SinterNormalSpecies::SinterNormalSpecies ( const SinterNormalSpecies p)

The default copy constructor.

Parameters
[in]thespecies that is copied

Definition at line 52 of file SinterNormalSpecies.cc.

References cohesionStiffness_, dissipation_, inverseSinterViscosity_, loadingStiffness_, penetrationDepthMax_, sinterAdhesion_, sinterRate_, sinterType_, and unloadingStiffnessMax_.

53 {
63 #ifdef DEBUG_CONSTRUCTOR
64  std::cout<<"SinterNormalSpecies::SinterNormalSpecies(const SinterNormalSpecies &p) finished"<<std::endl;
65 #endif
66 }
Mdouble unloadingStiffnessMax_
the maximum elastic constant (k_2^max) for plastic deformations
SINTERTYPE sinterType_
sinterType options determin how the rate of sintering, d(equilibriumOverlap)/dt is computed PARHAMI_M...
Mdouble loadingStiffness_
(normal) spring constant (k_1)
Mdouble cohesionStiffness_
the adhesive spring constant (k^c) for plastic deformations
Mdouble penetrationDepthMax_
the depth (relative to the normalized radius) at which k_2^max is used (phi_f)
Mdouble dissipation_
linear dissipation coefficient
SinterNormalSpecies::~SinterNormalSpecies ( )
virtual

The default destructor.

Definition at line 68 of file SinterNormalSpecies.cc.

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

Member Function Documentation

Mdouble SinterNormalSpecies::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.

Definition at line 221 of file SinterNormalSpecies.cc.

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

222 {
223  if (unloadingStiffnessMax_ / (.5 * mass) < mathsFunc::square(dissipation_ /mass)) {
224  std::cerr << "Dissipation too high; max. allowed " << sqrt(2.0 * unloadingStiffnessMax_ * mass) << std::endl;
225  exit(-1);
226  }
227  return 0.02 * constants::pi / std::sqrt(unloadingStiffnessMax_ / (.5 * mass) - mathsFunc::square(dissipation_ /mass));
228 }
Mdouble unloadingStiffnessMax_
the maximum elastic constant (k_2^max) for plastic deformations
T square(T val)
squares a number
Definition: ExtendedMath.h:91
const Mdouble pi
Definition: ExtendedMath.h:42
Mdouble dissipation_
linear dissipation coefficient
std::string SinterNormalSpecies::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")

Definition at line 113 of file SinterNormalSpecies.cc.

114 {
115  return "Sinter";
116 }
Mdouble SinterNormalSpecies::getCohesionStiffness ( ) const

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

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

Definition at line 172 of file SinterNormalSpecies.cc.

References cohesionStiffness_.

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

173 {
174  return cohesionStiffness_;
175 }
Mdouble cohesionStiffness_
the adhesive spring constant (k^c) for plastic deformations
Mdouble SinterNormalSpecies::getCollisionTime ( Mdouble  mass)

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

Definition at line 427 of file SinterNormalSpecies.cc.

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

428 {
429  if (mass <= 0)
430  {
431  std::cerr << "Error in getCollisionTime(" << mass << ") mass is not set or has an unexpected value, (getCollisionTime(" << mass << "))" << std::endl;
432  exit(-1);
433  }
434  if (loadingStiffness_ <= 0)
435  {
436  std::cerr << "Error in getCollisionTime(" << mass << ") stiffness=" << loadingStiffness_ << " is not set or has an unexpected value, (getCollisionTime(" << mass << "), with stiffness=" << loadingStiffness_ << ")" << std::endl;
437  exit(-1);
438  }
439  if (dissipation_ < 0)
440  {
441  std::cerr << "Error in getCollisionTime(" << mass << ") dissipation=" << dissipation_ << " is not set or has an unexpected value, (getCollisionTime(" << mass << "), with dissipation=" << dissipation_ << ")" << std::endl;
442  exit(-1);
443  }
444  Mdouble tosqrt = loadingStiffness_ / (.5 * mass) - mathsFunc::square(dissipation_ / mass);
445  if (tosqrt <= -1e-8*loadingStiffness_ / (.5 * mass))
446  {
447  std::cerr << "Error in getCollisionTime(" << mass << ") values for mass, stiffness and dissipation would lead to an overdamped system, (getCollisionTime(" << mass << "), with stiffness=" << loadingStiffness_ << " and dissipation=" << dissipation_ << ")" << std::endl;
448  exit(-1);
449  }
450  return constants::pi / std::sqrt(tosqrt);
451 }
Mdouble loadingStiffness_
(normal) spring constant (k_1)
double Mdouble
T square(T val)
squares a number
Definition: ExtendedMath.h:91
const Mdouble pi
Definition: ExtendedMath.h:42
Mdouble dissipation_
linear dissipation coefficient
Mdouble SinterNormalSpecies::getDissipation ( ) const

Allows the normal dissipation to be accessed.

Returns
the linear dissipation coefficient

Definition at line 343 of file SinterNormalSpecies.cc.

References dissipation_.

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

344 {
345  return dissipation_;
346 }
Mdouble dissipation_
linear dissipation coefficient
Mdouble SinterNormalSpecies::getInverseSinterViscosity ( ) const

Accesses inverseSinterViscosity_.

Returns
value of inverseSinterViscosity_

Definition at line 358 of file SinterNormalSpecies.cc.

References inverseSinterViscosity_.

Referenced by mix(), and setSinterForceAndTime().

359 {
361 }
Mdouble SinterNormalSpecies::getLoadingStiffness ( ) const

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

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

Definition at line 156 of file SinterNormalSpecies.cc.

References loadingStiffness_.

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

157 {
158  return loadingStiffness_;
159 }
Mdouble loadingStiffness_
(normal) spring constant (k_1)
Mdouble SinterNormalSpecies::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.

Definition at line 180 of file SinterNormalSpecies.cc.

References penetrationDepthMax_.

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

181 {
182  return penetrationDepthMax_;
183 }
Mdouble penetrationDepthMax_
the depth (relative to the normalized radius) at which k_2^max is used (phi_f)
Mdouble SinterNormalSpecies::getSinterAdhesion ( ) const

Accesses sinterAdhesion_.

Returns
value of sinterAdhesion_

Definition at line 351 of file SinterNormalSpecies.cc.

References sinterAdhesion_.

Referenced by SinterInteraction::computeNormalForce(), mix(), and setSinterForceAndTime().

352 {
353  return sinterAdhesion_;
354 }
Mdouble SinterNormalSpecies::getSinterRate ( ) const

Accesses sinterRate_.

Returns
value of sinterAdhesion_

Definition at line 365 of file SinterNormalSpecies.cc.

References sinterRate_.

Referenced by SinterInteraction::computeNormalForce().

366 {
367  return sinterRate_;
368 }
SINTERTYPE SinterNormalSpecies::getSinterType ( ) const

Definition at line 369 of file SinterNormalSpecies.cc.

References sinterType_.

Referenced by SinterInteraction::computeNormalForce().

370 {
371  return sinterType_;
372 }
SINTERTYPE sinterType_
sinterType options determin how the rate of sintering, d(equilibriumOverlap)/dt is computed PARHAMI_M...
std::function< double(double temperature)> SinterNormalSpecies::getTemperatureDependentSinterRate ( ) const

Definition at line 379 of file SinterNormalSpecies.cc.

References temperatureDependentSinterRate_.

Referenced by SinterInteraction::computeNormalForce().

380 {
382 }
std::function< double(double temperature)> temperatureDependentSinterRate_
double SinterNormalSpecies::getTemperatureDependentSinterRate ( double  temperature) const

Definition at line 374 of file SinterNormalSpecies.cc.

References temperatureDependentSinterRate_.

375 {
376  return temperatureDependentSinterRate_(temperature);
377 }
std::function< double(double temperature)> temperatureDependentSinterRate_
Mdouble SinterNormalSpecies::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.

Definition at line 164 of file SinterNormalSpecies.cc.

References unloadingStiffnessMax_.

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

165 {
166  return unloadingStiffnessMax_;
167 }
Mdouble unloadingStiffnessMax_
the maximum elastic constant (k_2^max) for plastic deformations
void SinterNormalSpecies::mix ( SinterNormalSpecies *const  S,
SinterNormalSpecies *const  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

Definition at line 123 of file SinterNormalSpecies.cc.

References BaseSpecies::average(), cohesionStiffness_, dissipation_, getCohesionStiffness(), getDissipation(), getInverseSinterViscosity(), getLoadingStiffness(), getPenetrationDepthMax(), getSinterAdhesion(), getUnloadingStiffnessMax(), inverseSinterViscosity_, loadingStiffness_, penetrationDepthMax_, sinterAdhesion_, and unloadingStiffnessMax_.

124 {
132 }
Mdouble unloadingStiffnessMax_
the maximum elastic constant (k_2^max) for plastic deformations
Mdouble loadingStiffness_
(normal) spring constant (k_1)
Mdouble cohesionStiffness_
the adhesive spring constant (k^c) for plastic deformations
Mdouble getSinterAdhesion() const
Accesses sinterAdhesion_.
Mdouble penetrationDepthMax_
the depth (relative to the normalized radius) at which k_2^max is used (phi_f)
Mdouble dissipation_
linear dissipation coefficient
Mdouble average(Mdouble a, Mdouble b)
defines the average of two variables by the harmonic mean.
Definition: BaseSpecies.cc:85
Mdouble getLoadingStiffness() const
Returns the loading stiffness of the linear plastic-viscoelastic normal force.
Mdouble getCohesionStiffness() const
Returns the cohesive stiffness of the linear plastic-viscoelastic normal force.
Mdouble getPenetrationDepthMax() const
Returns the maximum penetration depth of the linear plastic-viscoelastic normal force.
Mdouble getUnloadingStiffnessMax() const
Returns the maximum unloading stiffness of the linear plastic-viscoelastic normal force...
Mdouble getDissipation() const
Allows the normal dissipation to be accessed.
Mdouble getInverseSinterViscosity() const
Accesses inverseSinterViscosity_.
void SinterNormalSpecies::read ( std::istream &  is)
virtual

Reads the species properties from an input stream.

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

Implements BaseObject.

Definition at line 94 of file SinterNormalSpecies.cc.

References cohesionStiffness_, dissipation_, inverseSinterViscosity_, loadingStiffness_, penetrationDepthMax_, sinterAdhesion_, sinterRate_, sinterType_, and unloadingStiffnessMax_.

95 {
96  std::string dummy;
97  is >> dummy >> loadingStiffness_;
98  is >> dummy >> unloadingStiffnessMax_;
99  is >> dummy >> cohesionStiffness_;
100  is >> dummy >> penetrationDepthMax_;
101  is >> dummy >> dissipation_;
102  is >> dummy >> sinterAdhesion_;
103  is >> dummy >> inverseSinterViscosity_;
104  is >> dummy >> sinterRate_;
105  unsigned type;
106  is >> dummy >> type;
107  sinterType_ = (SINTERTYPE) type;
108 }
Mdouble unloadingStiffnessMax_
the maximum elastic constant (k_2^max) for plastic deformations
SINTERTYPE sinterType_
sinterType options determin how the rate of sintering, d(equilibriumOverlap)/dt is computed PARHAMI_M...
Mdouble loadingStiffness_
(normal) spring constant (k_1)
Mdouble cohesionStiffness_
the adhesive spring constant (k^c) for plastic deformations
Mdouble penetrationDepthMax_
the depth (relative to the normalized radius) at which k_2^max is used (phi_f)
Mdouble dissipation_
linear dissipation coefficient
void SinterNormalSpecies::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.

Definition at line 204 of file SinterNormalSpecies.cc.

References cohesionStiffness_.

Referenced by setPlasticParameters().

205 {
206  cohesionStiffness_ = cohesionStiffness;
207 }
Mdouble cohesionStiffness_
the adhesive spring constant (k^c) for plastic deformations
void SinterNormalSpecies::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?))

Definition at line 397 of file SinterNormalSpecies.cc.

398 {
399  if (eps==0.0) {
401  dissipation_ = std::sqrt(2.0 * mass * loadingStiffness_);
402  } else {
403  dissipation_ = -mass / tc * std::log(eps);
405  }
407 }
Mdouble unloadingStiffnessMax_
the maximum elastic constant (k_2^max) for plastic deformations
Mdouble loadingStiffness_
(normal) spring constant (k_1)
Mdouble log(Mdouble Power)
Definition: ExtendedMath.cc:97
T square(T val)
squares a number
Definition: ExtendedMath.h:91
const Mdouble pi
Definition: ExtendedMath.h:42
Mdouble dissipation_
linear dissipation coefficient
void SinterNormalSpecies::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.

Definition at line 234 of file SinterNormalSpecies.cc.

References dissipation_, ERROR, and logger.

Referenced by setLoadingStiffnessAndDissipation().

235 {
236  if (dissipation < 0)
237  {
238  logger(ERROR,"setDissipation(%)",dissipation);
239  exit(-1);
240  }
241  dissipation_ = dissipation;
242 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
Mdouble dissipation_
linear dissipation coefficient
void SinterNormalSpecies::setInverseSinterViscosity ( Mdouble  inverseSinterViscosity)

Sets inverseSinterViscosity_.

should be non-negative

Parameters
[in]thelinear dissipation coefficient of the linear plastic-viscoelastic normal force.

Definition at line 291 of file SinterNormalSpecies.cc.

References ERROR, inverseSinterViscosity_, logger, PARHAMI_MCKEEPING, and setSinterType().

Referenced by setSinterForceAndTime().

292 {
293  //assertOrDie(inverseSinterViscosity < 0, "inverseSinterViscosity should be non-negative!");
294  if (inverseSinterViscosity < 0)
295  {
296  logger(ERROR,"setInverseSinterViscosity(%)",inverseSinterViscosity);
297  }
299  inverseSinterViscosity_ = inverseSinterViscosity;
300 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
void setSinterType(SINTERTYPE sinterType)
Sets sinterRate_.
void SinterNormalSpecies::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.

Definition at line 188 of file SinterNormalSpecies.cc.

References loadingStiffness_.

Referenced by setLoadingStiffnessAndDissipation(), and setPlasticParameters().

189 {
190  loadingStiffness_ = loadingStiffness;
191 }
Mdouble loadingStiffness_
(normal) spring constant (k_1)
void SinterNormalSpecies::setLoadingStiffnessAndDissipation ( helpers::KAndDisp  new_)

Allows the spring and dissipation constants to be changed simultaneously.

Parameters
[in]new_a helper struct containing both the loading stiffness and the disssipation coefficient.

Definition at line 334 of file SinterNormalSpecies.cc.

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

335 {
336  setLoadingStiffness(new_.k);
337  setDissipation(new_.disp);
338 }
void setDissipation(Mdouble dissipation)
Sets the linear dissipation coefficient of the linear plastic-viscoelastic normal force...
Mdouble disp
Definition: Helpers.h:46
void setLoadingStiffness(Mdouble loadingStiffness)
Sets the loading stiffness of the linear plastic-viscoelastic normal force.
void SinterNormalSpecies::setParhamiMcKeeping ( Mdouble  alpha,
Mdouble  beta,
Mdouble  atomicVolume,
Mdouble  surfaceEnergy,
Mdouble  thicknessDiffusion,
Mdouble  activationEnergy,
Mdouble  temperature 
)

Sets the sinterAdhesion_ and inverseSinterViscosity_ based on the Parhami-McKeeping parameters.

should be non-negative

Parameters
[in]thelinear dissipation coefficient of the linear plastic-viscoelastic normal force.

Definition at line 317 of file SinterNormalSpecies.cc.

References mathsFunc::exp(), INFO, logger, PARHAMI_MCKEEPING, and constants::pi.

319 {
321  const Mdouble boltzmannConstant /*k_B*/ = 1.38064852e-23;
322  const Mdouble gasConstant /*R_g*/ = 8.314459848;
323  const Mdouble thicknessDiffusionVacancy /*DB*/ = thicknessDiffusion*exp(-activationEnergy/gasConstant/temperature);
324  const Mdouble diffusionParameter /*DeltaB*/ = atomicVolume/boltzmannConstant/temperature*thicknessDiffusionVacancy;
325  setInverseSinterViscosity(constants::pi / (2.0 * beta * diffusionParameter));
326  setSinterAdhesion(alpha/beta*constants::pi*surfaceEnergy);
327  logger(INFO,"Set sintering parameters: the adhesion force f_a=%*r, the rate of change of the plastic overlap is d(delta0)/dt=f_ep/(%*a^4)", getSinterAdhesion(),
329 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
Mdouble exp(Mdouble Exponent)
Definition: ExtendedMath.cc:78
Mdouble getSinterAdhesion() const
Accesses sinterAdhesion_.
double Mdouble
void setSinterAdhesion(Mdouble sinterAdhesion)
Sets sinterAdhesion_.
void setInverseSinterViscosity(Mdouble inverseSinterViscosity)
Sets inverseSinterViscosity_.
const Mdouble pi
Definition: ExtendedMath.h:42
void setSinterType(SINTERTYPE sinterType)
Sets sinterRate_.
Mdouble getInverseSinterViscosity() const
Accesses inverseSinterViscosity_.
void SinterNormalSpecies::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.

Definition at line 212 of file SinterNormalSpecies.cc.

References penetrationDepthMax_.

Referenced by setPlasticParameters().

213 {
214  penetrationDepthMax_ = penetrationDepthMax;
215 }
Mdouble penetrationDepthMax_
the depth (relative to the normalized radius) at which k_2^max is used (phi_f)
void SinterNormalSpecies::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.

Definition at line 140 of file SinterNormalSpecies.cc.

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

141 {
142  if (loadingStiffness <= 0 || unloadingStiffnessMax < loadingStiffness || cohesionStiffness < 0 || penetrationDepthMax < 0)
143  {
144  std::cerr << "Error: arguments of setPlasticParameters do not make sense" << std::endl;
145  exit(-1);
146  }
147  setLoadingStiffness(loadingStiffness);
148  setUnloadingStiffnessMax(unloadingStiffnessMax);
149  setCohesionStiffness(cohesionStiffness);
150  setPenetrationDepthMax(penetrationDepthMax);
151 }
void setCohesionStiffness(Mdouble cohesionStiffness)
Sets the cohesive stiffness of the linear plastic-viscoelastic normal force.
void setPenetrationDepthMax(Mdouble penetrationDepthMax)
Sets the maximum penetration depth of the linear plastic-viscoelastic normal force.
void setLoadingStiffness(Mdouble loadingStiffness)
Sets the loading stiffness of the linear plastic-viscoelastic normal force.
void setUnloadingStiffnessMax(Mdouble unloadingStiffnessMax)
Sets the maximum unloading stiffness of the linear plastic-viscoelastic normal force.
void SinterNormalSpecies::setSinterAdhesion ( Mdouble  sinterAdhesion)

Sets sinterAdhesion_.

should be non-negative

Parameters
[in]sinterAdhesion_

Definition at line 248 of file SinterNormalSpecies.cc.

References ERROR, logger, and sinterAdhesion_.

Referenced by setSinterForceAndTime().

249 {
250  if (sinterAdhesion < 0)
251  {
252  logger(ERROR,"setSinterAdhesion(%)",sinterAdhesion);
253  exit(-1);
254  }
255  sinterAdhesion_ = sinterAdhesion;
256 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
void SinterNormalSpecies::setSinterForceAndTime ( Mdouble  adhesionForce,
Mdouble  sinterTime,
Mdouble  radius 
)
Todo:
fix

Definition at line 302 of file SinterNormalSpecies.cc.

References getInverseSinterViscosity(), getSinterAdhesion(), INFO, logger, PARHAMI_MCKEEPING, setInverseSinterViscosity(), setSinterAdhesion(), and sinterType_.

303 {
306  setSinterAdhesion(adhesionForce/radius);
307  setInverseSinterViscosity(1.0/(93.75 * adhesionForce * sinterTime / std::pow(radius, 5)));
308  logger(INFO,"Set sintering parameters: adhesion force f_a=%*r, sinter viscosity is nu = contactRadius^4/%",
310 }
SINTERTYPE sinterType_
sinterType options determin how the rate of sintering, d(equilibriumOverlap)/dt is computed PARHAMI_M...
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
Mdouble getSinterAdhesion() const
Accesses sinterAdhesion_.
void setSinterAdhesion(Mdouble sinterAdhesion)
Sets sinterAdhesion_.
void setInverseSinterViscosity(Mdouble inverseSinterViscosity)
Sets inverseSinterViscosity_.
Mdouble getInverseSinterViscosity() const
Accesses inverseSinterViscosity_.
void SinterNormalSpecies::setSinterRate ( Mdouble  sinterRate)

Sets sinterRate_.

should be non-negative

Parameters
[in]sinterAdhesion_

Definition at line 261 of file SinterNormalSpecies.cc.

References ERROR, logger, and sinterRate_.

262 {
263  if (sinterRate < 0)
264  {
265  logger(ERROR,"setSinterRate(%)",sinterRate);
266  }
267  sinterRate_ = sinterRate;
268 }
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
void SinterNormalSpecies::setSinterType ( SINTERTYPE  sinterType)

Sets sinterRate_.

should be non-negative

Parameters
[in]sinterAdhesion_

Definition at line 274 of file SinterNormalSpecies.cc.

References CONSTANT_RATE, ERROR, INFO, logger, PARHAMI_MCKEEPING, sinterType_, and TEMPERATURE_DEPENDENT_FRENKEL.

Referenced by setInverseSinterViscosity().

275 {
276  sinterType_ = sinterType;
278  logger(INFO, "Sintertype set to CONSTANT_RATE");
280  logger(INFO, "Sintertype set to PARHAMI_MCKEEPING");
282  logger(INFO, "Sintertype set to TEMPERATURE_DEPENDENT_FRENKEL");
283  else
284  logger(ERROR, "Sintertype not understood");
285 }
SINTERTYPE sinterType_
sinterType options determin how the rate of sintering, d(equilibriumOverlap)/dt is computed PARHAMI_M...
Logger< MERCURY_LOGLEVEL > logger("MercuryKernel")
void SinterNormalSpecies::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}$

Definition at line 415 of file SinterNormalSpecies.cc.

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

416 {
417  loadingStiffness_ = stiffness;
418  if (eps==0.0) {
419  dissipation_ = std::sqrt(2.0 * mass * stiffness);
420  } else {
421  dissipation_ = -std::sqrt(2.0 * mass * stiffness / (constants::sqr_pi + mathsFunc::square(log(eps)))) * log(eps);
422  }
424 }
Mdouble unloadingStiffnessMax_
the maximum elastic constant (k_2^max) for plastic deformations
Mdouble loadingStiffness_
(normal) spring constant (k_1)
Mdouble log(Mdouble Power)
Definition: ExtendedMath.cc:97
T square(T val)
squares a number
Definition: ExtendedMath.h:91
Mdouble dissipation_
linear dissipation coefficient
const Mdouble sqr_pi
Definition: ExtendedMath.h:44
void SinterNormalSpecies::setTemperatureDependentSinterRate ( std::function< double(double temperature)>  temperatureDependentSinterRate)

Definition at line 384 of file SinterNormalSpecies.cc.

References temperatureDependentSinterRate_.

385 {
386  temperatureDependentSinterRate_ = temperatureDependentSinterRate;
387 }
std::function< double(double temperature)> temperatureDependentSinterRate_
void SinterNormalSpecies::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.

Definition at line 196 of file SinterNormalSpecies.cc.

References unloadingStiffnessMax_.

Referenced by setPlasticParameters().

197 {
198  unloadingStiffnessMax_ = unloadingStiffnessMax;
199 }
Mdouble unloadingStiffnessMax_
the maximum elastic constant (k_2^max) for plastic deformations
void SinterNormalSpecies::write ( std::ostream &  os) const
virtual

Writes the species properties to an output stream.

Parameters
[out]outputstream (typically the restart file)

Implements BaseObject.

Definition at line 78 of file SinterNormalSpecies.cc.

References cohesionStiffness_, dissipation_, inverseSinterViscosity_, loadingStiffness_, penetrationDepthMax_, sinterAdhesion_, sinterRate_, sinterType_, and unloadingStiffnessMax_.

79 {
80  os << " loadingStiffness " << loadingStiffness_;
81  os << " maxUnloadingStiffness " << unloadingStiffnessMax_;
82  os << " cohesionStiffness " << cohesionStiffness_;
83  os << " maxPenetration " << penetrationDepthMax_;
84  os << " dissipation " << dissipation_;
85  os << " sinterAdhesion " << sinterAdhesion_;
86  os << " inverseSinterViscosity " << inverseSinterViscosity_;
87  os << " sinterRate_ " << sinterRate_;
88  os << " sinterType " << (unsigned) sinterType_;
89 }
Mdouble unloadingStiffnessMax_
the maximum elastic constant (k_2^max) for plastic deformations
SINTERTYPE sinterType_
sinterType options determin how the rate of sintering, d(equilibriumOverlap)/dt is computed PARHAMI_M...
Mdouble loadingStiffness_
(normal) spring constant (k_1)
Mdouble cohesionStiffness_
the adhesive spring constant (k^c) for plastic deformations
Mdouble penetrationDepthMax_
the depth (relative to the normalized radius) at which k_2^max is used (phi_f)
Mdouble dissipation_
linear dissipation coefficient

Member Data Documentation

Mdouble SinterNormalSpecies::cohesionStiffness_
private

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

Definition at line 212 of file SinterNormalSpecies.h.

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

Mdouble SinterNormalSpecies::dissipation_
private
Mdouble SinterNormalSpecies::inverseSinterViscosity_
private

Determines sinter viscosity nu = contactRadius^4/inverseSinterViscosity_ in sinter rate: d(delta0)/dt = (fa+fep)/nu

Definition at line 228 of file SinterNormalSpecies.h.

Referenced by getInverseSinterViscosity(), mix(), read(), setInverseSinterViscosity(), SinterNormalSpecies(), and write().

Mdouble SinterNormalSpecies::loadingStiffness_
private
Mdouble SinterNormalSpecies::penetrationDepthMax_
private

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

Definition at line 215 of file SinterNormalSpecies.h.

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

Mdouble SinterNormalSpecies::sinterAdhesion_
private

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

Definition at line 223 of file SinterNormalSpecies.h.

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

Mdouble SinterNormalSpecies::sinterRate_
private

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

Definition at line 233 of file SinterNormalSpecies.h.

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

SINTERTYPE SinterNormalSpecies::sinterType_
private

sinterType options determin how the rate of sintering, d(equilibriumOverlap)/dt is computed PARHAMI_MCKEEPING: sinter rate given by (sinterAdhesion+normalForce)/sinterViscosity CONSTANT_RATE: sinter rate given by sinterRate

Definition at line 240 of file SinterNormalSpecies.h.

Referenced by getSinterType(), read(), setSinterForceAndTime(), setSinterType(), SinterNormalSpecies(), and write().

std::function<double(double temperature)> SinterNormalSpecies::temperatureDependentSinterRate_ = [this] (double temperature) {return sinterRate_;}
private
Mdouble SinterNormalSpecies::unloadingStiffnessMax_
private

the maximum elastic constant (k_2^max) for plastic deformations

Definition at line 209 of file SinterNormalSpecies.h.

Referenced by computeTimeStep(), getUnloadingStiffnessMax(), mix(), read(), setStiffnessAndRestitutionCoefficient(), setUnloadingStiffnessMax(), SinterNormalSpecies(), and write().


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