HG-MD  1
Public Member Functions | Public Attributes | Friends
CSpecies Class Reference

Stores properties of the particles and the contact models such as the elastic modulus. More...

#include <CSpecies.h>

List of all members.

Public Member Functions

 CSpecies ()
 CSpecies (const CSpecies &S)
void set_k (Mdouble new_k)
 Allows the spring constant to be changed.
Mdouble get_k () const
 Allows the spring constant to be accessed.
void set_kt (Mdouble new_kt)
 Allows the spring constant to be changed.
Mdouble get_kt () const
 Allows the spring constant to be accessed.
void set_krolling (Mdouble new_k)
 Allows the spring constant to be changed.
Mdouble get_krolling () const
 Allows the spring constant to be accessed.
void set_ktorsion (Mdouble new_k)
 Allows the spring constant to be changed.
Mdouble get_ktorsion () const
 Allows the spring constant to be accessed.
void set_rho (Mdouble new_rho)
Mdouble get_rho () const
 Allows the density to be accessed.
void set_dispt (Mdouble new_dispt)
 Allows the tangential viscosity to be changed.
Mdouble get_dispt () const
 Allows the tangential viscosity to be accessed.
void set_disprolling (Mdouble new_disprolling)
 Allows the tangential viscosity to be changed.
Mdouble get_disprolling () const
 Allows the tangential viscosity to be accessed.
void set_disptorsion (Mdouble new_disptorsion)
 Allows the tangential viscosity to be changed.
Mdouble get_disptorsion () const
 Allows the tangential viscosity to be accessed.
void set_dissipation (Mdouble new_disp)
 Allows the normal dissipation to be changed.
Mdouble get_dissipation () const
 Allows the normal dissipation to be accessed.
void set_mu (Mdouble new_mu)
 Allows the (dynamic) Coulomb friction coefficient to be changed; also sets mu_s by default.
Mdouble get_mu () const
 Allows the (dynamic) Coulomb friction coefficient to be accessed.
void set_mus (Mdouble new_mu)
 Allows the static Coulomb friction coefficient to be changed.
Mdouble get_mus () const
 Allows the static Coulomb friction coefficient to be accessed.
void set_murolling (Mdouble new_murolling)
 Allows the (dynamic) Coulomb friction coefficient to be changed; also sets murolling_s by default.
Mdouble get_murolling () const
 Allows the (dynamic) Coulomb friction coefficient to be accessed.
void set_musrolling (Mdouble new_mu)
 Allows the static Coulomb friction coefficient to be changed.
Mdouble get_musrolling () const
 Allows the static Coulomb friction coefficient to be accessed.
void set_mutorsion (Mdouble new_mu)
 Allows the (dynamic) Coulomb friction coefficient to be changed; also sets mu_s by default.
Mdouble get_mutorsion () const
 Allows the (dynamic) Coulomb friction coefficient to be accessed.
void set_mustorsion (Mdouble new_mu)
 Allows the static Coulomb friction coefficient to be changed.
Mdouble get_mustorsion () const
 Allows the static Coulomb friction coefficient to be accessed.
void set_dim_particle (int new_dim)
 Allows the dimension of the particle (f.e. for mass) to be changed.
int get_dim_particle () const
 Allows the dimension of the particle (f.e. for mass) to be accessed.
void read (std::istream &is)
void print (std::ostream &os)
 Outputs species.
Mdouble get_collision_time (Mdouble mass)
 Calculates collision time for two copies of a particle of given disp, k, mass.
Mdouble get_restitution_coefficient (Mdouble mass)
 Calculates restitution coefficient for two copies of given disp, k, mass.
Mdouble get_maximum_velocity (Mdouble radius, Mdouble mass)
 Calculates the maximum velocity allowed for a collision of two copies of P (for higher velocities particles could pass through each other)
void set_k_and_restitution_coefficient (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.
void set_collision_time_and_restitution_coefficient (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.
void set_collision_time_and_restitution_coefficient (Mdouble tc, Mdouble eps, Mdouble mass1, Mdouble mass2)
 Set k, disp such that is matches a given tc and eps for a collision of two different masses.
void set_collision_time_and_normal_and_tangential_restitution_coefficient (Mdouble tc, Mdouble eps, Mdouble beta, Mdouble mass)
 Sets k, disp, kt, dispt such that it matches a given tc and eps for a collision of two particles of masses m0,m1.
void set_collision_time_and_normal_and_tangential_restitution_coefficient_nodispt (Mdouble tc, Mdouble eps, Mdouble beta, Mdouble mass)
 Sets k, disp, kt, dispt such that it matches a given tc and eps for a collision of two particles of masses m0,m1.
void set_collision_time_and_normal_and_tangential_restitution_coefficient (Mdouble tc, Mdouble eps, Mdouble beta, Mdouble mass1, Mdouble mass2)
 Sets k, disp, kt, dispt such that it matches a given tc and eps for a collision of two particles of equal mass m.
void set_collision_time_and_normal_and_tangential_restitution_coefficient_nodispt (Mdouble tc, Mdouble eps, Mdouble beta, Mdouble mass1, Mdouble mass2)
 Sets k, disp, kt, dispt such that it matches a given tc and eps for a collision of two particles of equal mass m.
Mdouble get_average (Mdouble a, Mdouble b)
void mix (CSpecies &S0, CSpecies &S1)
 create values for mixed species
void set_plastic_k1_k2max_kc_depth (Mdouble k1_, Mdouble k2max_, Mdouble kc_, Mdouble depth_)
 Acccess functions for the plastic model.
Mdouble get_k1 () const
Mdouble get_k2max () const
Mdouble get_kc () const
Mdouble get_depth () const
void set_k1 (Mdouble new_)
void set_k2max (Mdouble new_)
void set_kc (Mdouble new_)
void set_depth (Mdouble new_)
Mdouble get_plastic_dt (Mdouble mass)
 Calculates collision time for stiffest spring constant, divides by 50.

Public Attributes

Mdouble k
Mdouble kt
Mdouble krolling
Mdouble ktorsion
Mdouble disp
Mdouble dispt
Mdouble disprolling
Mdouble disptorsion
Mdouble mu
Mdouble mus
Mdouble murolling
Mdouble musrolling
Mdouble mutorsion
Mdouble mustorsion
Mdouble rho
int dim_particle
vector< CSpeciesMixedSpecies
Mdouble k2max
Mdouble kc
Mdouble depth

Friends

std::ostream & operator<< (std::ostream &os, const CSpecies &s)
 Writes all species data.
std::istream & operator>> (std::istream &is, CSpecies &s)
 Reads all species data.

Detailed Description

Stores properties of the particles and the contact models such as the elastic modulus.

As particle species are distinguished by their species Particle::indSpecies. Based on this index, different particle and contact properties apply which are stored in MD::Species. The contact model for interspecies interactions are defined in CSpecies::MixedSpecies

Bug:
Segfault, caused by: In MD::add_Species(...) (MD.cc, line 1298), if the line Species.push_back(S); causes a reallocation of internal storage, then the existing elements of Species are copy-constructed into the new memory. However, the copy constructor of CSpecies (CSpecies.h, line 35), does not copy the MixedSpecies vector, so if we call add_Species more than once (meaning that for the second and subsequent calls, existing elements in the Species vector have non-empty MixedSpecies vectors) we are liable to lose this data. Fix: Copy the MixedSpecies member in the copy constructor (CSpecies.h, line 40-ish): MixedSpecies = S.MixedSpecies; However, it looks to me as though if this line is added, the copy constructor could be replaced by the default bitwise copy constructor, so perhaps the copying of MixedSpecies has been left out for a reason...

Constructor & Destructor Documentation

CSpecies::CSpecies ( ) [inline]

References depth, dim_particle, disp, disprolling, dispt, disptorsion, k, k2max, kc, krolling, kt, ktorsion, mu, murolling, mus, musrolling, mustorsion, mutorsion, and rho.

                   {
                k = 0; 
                kt = 0;
                krolling = 0;
                ktorsion = 0;
                disp = 0;
                dispt = 0;
                disprolling = 0;
                disptorsion = 0;
                mu = 0;
                mus = 0;
                mu = 0;
                murolling = 0;
                musrolling = 0;
                mutorsion = 0;
                mustorsion = 0;
                k2max = 0;
                kc = 0;
                depth = 0;
                rho = 0;
                dim_particle = 0;
        }
CSpecies::CSpecies ( const CSpecies S) [inline]

Member Function Documentation

Mdouble CSpecies::get_average ( Mdouble  a,
Mdouble  b 
) [inline]

Referenced by mix().

                                                  {
                return (a+b) ? (2.*(a*b)/(a+b)) : 0;
        }

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

References disp, k, constants::pi, and sqr.

Referenced by get_restitution_coefficient().

                                                {
                if(mass<=0)     {cerr << "Error in get_collision_time(Mdouble mass) mass is not set or has an unexpected value, (get_collision_time("<<mass<<"))"<< endl; exit(-1);}
                if(k<=0)                {cerr << "Error in get_collision_time(Mdouble mass) stiffness is not set or has an unexpected value, (get_collision_time("<<mass<<"), with stiffness="<<k<<")"<< endl; exit(-1);}
                if(disp<0)              {cerr << "Error in get_collision_time(Mdouble mass) dissipation is not set or has an unexpected value, (get_collision_time("<<mass<<"), with dissipation="<<disp<<")"<< endl; exit(-1);}
                Mdouble tosqrt=k/(.5*mass) - sqr(disp/mass);
                //~ cout<<"tosqrt "<<tosqrt<<" 1e "<<k/(.5*mass)<<" 2e "<<sqr(disp/mass)<<endl;
                if (tosqrt<=0)  {cerr << "Error in get_collision_time(Mdouble mass) values for mass, stiffness and dissipation would leads to an over damped system, (get_collision_time("<<mass<<"), with stiffness="<<k<<" and dissipation="<<disp<<")"<< endl; exit(-1);}
                return constants::pi / sqrt( tosqrt );
        }
Mdouble CSpecies::get_depth ( ) const [inline]

References depth.

Referenced by CSpecies(), and mix().

{return depth;}
int CSpecies::get_dim_particle ( ) const [inline]

Allows the dimension of the particle (f.e. for mass) to be accessed.

References dim_particle.

Referenced by CSpecies().

{return dim_particle;}
Mdouble CSpecies::get_disprolling ( ) const [inline]

Allows the tangential viscosity to be accessed.

References disprolling.

Referenced by CSpecies(), and mix().

{return disprolling;}
Mdouble CSpecies::get_dispt ( ) const [inline]

Allows the tangential viscosity to be accessed.

References dispt.

Referenced by CSpecies(), and mix().

{return dispt;}
Mdouble CSpecies::get_disptorsion ( ) const [inline]

Allows the tangential viscosity to be accessed.

References disptorsion.

Referenced by CSpecies(), and mix().

{return disptorsion;}
Mdouble CSpecies::get_dissipation ( ) const [inline]

Allows the normal dissipation to be accessed.

References disp.

Referenced by CSpecies(), and mix().

{return disp;}
Mdouble CSpecies::get_k ( ) const [inline]
Mdouble CSpecies::get_k1 ( ) const [inline]

References k.

{return k;}
Mdouble CSpecies::get_k2max ( ) const [inline]

References k2max.

Referenced by CSpecies(), and mix().

{return k2max;}
Mdouble CSpecies::get_kc ( ) const [inline]

References kc.

Referenced by CSpecies(), and mix().

{return kc;}
Mdouble CSpecies::get_krolling ( ) const [inline]

Allows the spring constant to be accessed.

References krolling.

Referenced by CSpecies(), and mix().

{return krolling;}
Mdouble CSpecies::get_kt ( ) const [inline]

Allows the spring constant to be accessed.

References kt.

Referenced by CSpecies(), mix(), and set_collision_time_and_normal_and_tangential_restitution_coefficient().

{return kt;}
Mdouble CSpecies::get_ktorsion ( ) const [inline]

Allows the spring constant to be accessed.

References ktorsion.

Referenced by CSpecies(), and mix().

{return ktorsion;}
Mdouble CSpecies::get_maximum_velocity ( Mdouble  radius,
Mdouble  mass 
) [inline]

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

References k.

{return radius * sqrt(k/(.5*mass));}
Mdouble CSpecies::get_mu ( ) const [inline]

Allows the (dynamic) Coulomb friction coefficient to be accessed.

References mu.

Referenced by CSpecies(), and mix().

{return mu;}
Mdouble CSpecies::get_murolling ( ) const [inline]

Allows the (dynamic) Coulomb friction coefficient to be accessed.

References murolling.

Referenced by CSpecies(), and mix().

{return murolling;}
Mdouble CSpecies::get_mus ( ) const [inline]

Allows the static Coulomb friction coefficient to be accessed.

References mus.

Referenced by CSpecies(), and mix().

{return mus;}
Mdouble CSpecies::get_musrolling ( ) const [inline]

Allows the static Coulomb friction coefficient to be accessed.

References musrolling.

Referenced by CSpecies(), and mix().

{return musrolling;}
Mdouble CSpecies::get_mustorsion ( ) const [inline]

Allows the static Coulomb friction coefficient to be accessed.

References mustorsion.

Referenced by CSpecies(), and mix().

{return mustorsion;}
Mdouble CSpecies::get_mutorsion ( ) const [inline]

Allows the (dynamic) Coulomb friction coefficient to be accessed.

References mutorsion.

Referenced by CSpecies(), and mix().

{return mutorsion;}

Calculates collision time for stiffest spring constant, divides by 50.

References disp, k2max, constants::pi, and sqr.

                                            {
                return 0.02 * constants::pi / sqrt( k2max/(.5*mass) - sqr(disp/mass) );
        }

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

References disp, and get_collision_time().

{return exp(-disp/mass*get_collision_time(mass));}
Mdouble CSpecies::get_rho ( ) const [inline]

Allows the density to be accessed.

References rho.

Referenced by CSpecies().

{return rho;}
void CSpecies::mix ( CSpecies S0,
CSpecies S1 
) [inline]

create values for mixed species

References depth, dim_particle, disp, disprolling, dispt, disptorsion, get_average(), get_depth(), get_disprolling(), get_dispt(), get_disptorsion(), get_dissipation(), get_k(), get_k2max(), get_kc(), get_krolling(), get_kt(), get_ktorsion(), get_mu(), get_murolling(), get_mus(), get_musrolling(), get_mustorsion(), get_mutorsion(), k, k2max, kc, krolling, kt, ktorsion, mu, murolling, mus, musrolling, mustorsion, mutorsion, and rho.

void CSpecies::print ( std::ostream &  os) [inline]

Outputs species.

References depth, dim_particle, disp, disprolling, dispt, disptorsion, k, k2max, kc, krolling, kt, ktorsion, mu, murolling, mus, musrolling, mustorsion, mutorsion, and rho.

                                   {
                os      << "k:" << k 
                        << ", disp:" << disp 
                        << ", kt:" << kt;
                if (krolling) os << ", krolling: " << krolling;
                if (ktorsion) os << ", ktorsion: " << ktorsion;
                os  << ", dispt: " << dispt;
                if (disprolling) os << ", disprolling: " << disprolling;
                if (disptorsion) os << ", disptorsion: " << disptorsion;
                os << ", mu: " << mu;
                if (mu!=mus) os << ", mus: " << mus;
                if (murolling) os << ", murolling: " << murolling;
                if (murolling!=musrolling) os << ", musrolling: " << musrolling;
                if (mutorsion) os << ", mutorsion: " << mutorsion;
                if (mutorsion!=mustorsion) os << ", mustorsion: " << mustorsion;
                os      << ", k2max:" << k2max 
                        << ", kc:" << kc 
                        << ", depth:" << depth;
                if (dim_particle) {
                        os      << ", rho:" << rho 
                                << ", dim_particle:" << dim_particle;
                } else {
                        os      << " (mixed)";
                }
        }
void CSpecies::read ( std::istream &  is) [inline]

References depth, dim_particle, disp, disprolling, dispt, disptorsion, k, k2max, kc, krolling, kt, ktorsion, mu, murolling, mus, musrolling, mustorsion, mutorsion, and rho.

                                   { 
                string dummy;
                is      >> dummy >> k 
                        >> dummy >> disp 
                        >> dummy >> kt >> dummy;
                if (!strcmp(dummy.c_str(),"krolling")) is >> krolling >> dummy; else krolling = 0;
                if (!strcmp(dummy.c_str(),"ktorsion")) is >> ktorsion >> dummy; else ktorsion = 0;
                is      >> dispt >> dummy;
                if (!strcmp(dummy.c_str(),"disprolling")) is >> disprolling >> dummy; else disprolling = 0;
                if (!strcmp(dummy.c_str(),"disptorsion")) is >> disptorsion >> dummy; else disptorsion = 0;
                is      >> mu >> dummy;
                if (!strcmp(dummy.c_str(),"mus")) is >> mus >> dummy; else mus=mu;
                if (!strcmp(dummy.c_str(),"murolling"))  is >> murolling  >> dummy; else murolling  = 0;
                if (!strcmp(dummy.c_str(),"musrolling")) is >> musrolling >> dummy; else musrolling = murolling;
                if (!strcmp(dummy.c_str(),"mutorsion"))  is >> mutorsion  >> dummy; else mutorsion  = 0;
                if (!strcmp(dummy.c_str(),"mustorsion")) is >> mustorsion >> dummy; else mustorsion = mutorsion;

                // checks if plastic model parameters are included (for backward compability)
                if (!strcmp(dummy.c_str(),"k2max")) {
                        is >> k2max 
                                >> dummy >> kc 
                                >> dummy >> depth
                                >> dummy;
                }
                is      >> rho 
                        >> dummy >> dim_particle;
                //read rest of line
                getline(is,dummy);
        }

Sets k, disp, kt, dispt such that it matches a given tc and eps for a collision of two particles of masses m0,m1.

References get_k(), get_kt(), constants::pi, set_collision_time_and_restitution_coefficient(), set_dispt(), set_kt(), and sqr.

Referenced by set_collision_time_and_normal_and_tangential_restitution_coefficient().

                                                                                                                                      {
                set_collision_time_and_restitution_coefficient(tc,eps,mass);
                //from Deen...Kuipers2006, eq. 43 and 30
                set_kt(2.0/7.0*get_k()*(sqr(constants::pi)+sqr(log(beta)))/(sqr(constants::pi)+sqr(log(eps))));
                if (beta) set_dispt(-2*log(beta)*sqrt(1.0/7.0*mass*get_kt()/(sqr(constants::pi)+sqr(log(beta)))));
                else set_dispt(2.*sqrt(1.0/7.0*mass*get_kt()));
        }

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

References set_collision_time_and_normal_and_tangential_restitution_coefficient().

                {
                        Mdouble reduced_mass = mass1*mass2/(mass1+mass2);
                        set_collision_time_and_normal_and_tangential_restitution_coefficient(tc,eps,beta,2.0*reduced_mass);
                }       

Sets k, disp, kt, dispt such that it matches a given tc and eps for a collision of two particles of masses m0,m1.

References get_k(), constants::pi, set_collision_time_and_restitution_coefficient(), set_dispt(), set_kt(), and sqr.

Referenced by set_collision_time_and_normal_and_tangential_restitution_coefficient_nodispt().

                                                                                                                                              {
                set_collision_time_and_restitution_coefficient(tc,eps,mass);
                //from BeckerSchwagerPoeschel2008, eq. 56
                set_kt(2.0/7.0*get_k()*sqr(acos(-beta)/constants::pi));
                set_dispt(0);
        }

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

References set_collision_time_and_normal_and_tangential_restitution_coefficient_nodispt().

                {
                        Mdouble reduced_mass = mass1*mass2/(mass1+mass2);
                        set_collision_time_and_normal_and_tangential_restitution_coefficient_nodispt(tc,eps,beta,2.0*reduced_mass);
                }       

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

References disp, k, constants::pi, and sqr.

Referenced by set_collision_time_and_normal_and_tangential_restitution_coefficient(), set_collision_time_and_normal_and_tangential_restitution_coefficient_nodispt(), and set_collision_time_and_restitution_coefficient().

                                                                                                  {
                disp = - mass / tc * log(eps);
                k = .5 * mass * (sqr(constants::pi/tc) + sqr(disp/mass));
        }
void CSpecies::set_collision_time_and_restitution_coefficient ( Mdouble  tc,
Mdouble  eps,
Mdouble  mass1,
Mdouble  mass2 
) [inline]

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 set_collision_time_and_restitution_coefficient(Mdouble tc, Mdouble eps, Mdouble mass)

References set_collision_time_and_restitution_coefficient().

                {
                        Mdouble reduced_mass = mass1*mass2/(mass1+mass2);
                        set_collision_time_and_restitution_coefficient(tc,eps,2.0*reduced_mass);
                }       
void CSpecies::set_depth ( Mdouble  new_) [inline]

References depth.

Referenced by set_plastic_k1_k2max_kc_depth().

{depth = new_;}
void CSpecies::set_dim_particle ( int  new_dim) [inline]

Allows the dimension of the particle (f.e. for mass) to be changed.

References dim_particle.

{if (new_dim>=1 && new_dim<=3) dim_particle = new_dim; else { cerr << "Error in set_dim_particle" << endl; exit(-1); }}
void CSpecies::set_disprolling ( Mdouble  new_disprolling) [inline]

Allows the tangential viscosity to be changed.

References disprolling.

{if (new_disprolling>=0) disprolling = new_disprolling; else { cerr << "Error in set_disprolling" << endl; exit(-1); }}
void CSpecies::set_dispt ( Mdouble  new_dispt) [inline]

Allows the tangential viscosity to be changed.

References dispt.

Referenced by set_collision_time_and_normal_and_tangential_restitution_coefficient(), and set_collision_time_and_normal_and_tangential_restitution_coefficient_nodispt().

{if (new_dispt>=0) dispt = new_dispt; else { cerr << "Error in set_dispt" << endl; exit(-1); }}
void CSpecies::set_disptorsion ( Mdouble  new_disptorsion) [inline]

Allows the tangential viscosity to be changed.

References disptorsion.

{if (new_disptorsion>=0) disptorsion = new_disptorsion; else { cerr << "Error in set_disptorsion" << endl; exit(-1); }}
void CSpecies::set_dissipation ( Mdouble  new_disp) [inline]

Allows the normal dissipation to be changed.

References disp.

{if(new_disp>=0) disp=new_disp; else { cerr << "Error in set_dissipation(" << new_disp << ")" << endl; exit(-1); }}
void CSpecies::set_k ( Mdouble  new_k) [inline]

Allows the spring constant to be changed.

References k.

{if(new_k>=0) k=new_k; else { cerr << "Error in set_k" << endl; exit(-1); }}
void CSpecies::set_k1 ( Mdouble  new_) [inline]

References k.

Referenced by set_plastic_k1_k2max_kc_depth().

{k = new_;}
void CSpecies::set_k2max ( Mdouble  new_) [inline]

References k2max.

Referenced by set_plastic_k1_k2max_kc_depth().

{k2max = new_;}
void CSpecies::set_k_and_restitution_coefficient ( Mdouble  k_,
Mdouble  eps,
Mdouble  mass 
) [inline]

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

References disp, k, sqr, and constants::sqrt_pi.

                                                                                     {
                k = k_;
                disp = - sqrt(2.0*mass*k/(constants::sqrt_pi+sqr(log(eps)))) * log(eps);
        }
void CSpecies::set_kc ( Mdouble  new_) [inline]

References kc.

Referenced by set_plastic_k1_k2max_kc_depth().

{kc = new_;}
void CSpecies::set_krolling ( Mdouble  new_k) [inline]

Allows the spring constant to be changed.

References krolling.

{if(new_k>=0) {krolling=new_k;} else { cerr << "Error in set_krolling" << endl; exit(-1); }}
void CSpecies::set_kt ( Mdouble  new_kt) [inline]

Allows the spring constant to be changed.

References kt.

Referenced by set_collision_time_and_normal_and_tangential_restitution_coefficient(), and set_collision_time_and_normal_and_tangential_restitution_coefficient_nodispt().

{if(new_kt>=0) {kt=new_kt;} else { cerr << "Error in set_kt" << endl; exit(-1); }}
void CSpecies::set_ktorsion ( Mdouble  new_k) [inline]

Allows the spring constant to be changed.

References ktorsion.

{if(new_k>=0) {ktorsion=new_k;} else { cerr << "Error in set_ktorsion" << endl; exit(-1); }}
void CSpecies::set_mu ( Mdouble  new_mu) [inline]

Allows the (dynamic) Coulomb friction coefficient to be changed; also sets mu_s by default.

References mu, and mus.

{if (new_mu>=0) {mu = new_mu; mus = mu;} else { cerr << "Error in set_mu" << endl; exit(-1); }}
void CSpecies::set_murolling ( Mdouble  new_murolling) [inline]

Allows the (dynamic) Coulomb friction coefficient to be changed; also sets murolling_s by default.

References murolling, and musrolling.

{if (new_murolling>=0) {murolling = new_murolling; musrolling = murolling;} else { cerr << "Error in set_murolling" << endl; exit(-1); }}
void CSpecies::set_mus ( Mdouble  new_mu) [inline]

Allows the static Coulomb friction coefficient to be changed.

References mus.

{if (new_mu>=0) {mus = new_mu;} else { cerr << "Error in set_mus" << endl; exit(-1); }}
void CSpecies::set_musrolling ( Mdouble  new_mu) [inline]

Allows the static Coulomb friction coefficient to be changed.

References musrolling.

{if (new_mu>=0) {musrolling = new_mu;} else { cerr << "Error in set_musrolling" << endl; exit(-1); }}
void CSpecies::set_mustorsion ( Mdouble  new_mu) [inline]

Allows the static Coulomb friction coefficient to be changed.

References mustorsion.

{if (new_mu>=0) {mustorsion = new_mu;} else { cerr << "Error in set_mustorsion" << endl; exit(-1); }}
void CSpecies::set_mutorsion ( Mdouble  new_mu) [inline]

Allows the (dynamic) Coulomb friction coefficient to be changed; also sets mu_s by default.

References mustorsion, and mutorsion.

{if (new_mu>=0) {mutorsion = new_mu; mustorsion = mutorsion;} else { cerr << "Error in set_mutorsion" << endl; exit(-1); }}
void CSpecies::set_plastic_k1_k2max_kc_depth ( Mdouble  k1_,
Mdouble  k2max_,
Mdouble  kc_,
Mdouble  depth_ 
) [inline]

Acccess functions for the plastic model.

References set_depth(), set_k1(), set_k2max(), and set_kc().

                                                                                                    {
                set_k1(k1_);
                set_k2max(k2max_);
                set_kc(kc_);
                set_depth(depth_);
        }       
void CSpecies::set_rho ( Mdouble  new_rho) [inline]
Todo:
recalculate masses when setting dim_particle or rho Allows the density to be changed

References rho.

{if(new_rho>=0) rho=new_rho; else { cerr << "Error in set_rho" << endl; exit(-1); }}

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  os,
const CSpecies s 
) [friend]

Writes all species data.

                                                                                {
                os << "k " << s.k 
                        << " disp " << s.disp 
                        << " kt " << s.kt;
                if (s.krolling) os << " krolling " << s.krolling;
                if (s.ktorsion) os << " ktorsion " << s.ktorsion;
                os  << " dispt " << s.dispt;
                if (s.disprolling) os << " disprolling " << s.disprolling;
                if (s.disptorsion) os << " disptorsion " << s.disptorsion;
                os << " mu " << s.mu;
                //optional output
                if (s.mu!=s.mus) os << " mus " << s.mus;
                if (s.murolling) os << " murolling " << s.murolling;
                if (s.murolling!=s.musrolling) os << " musrolling " << s.musrolling;
                if (s.mutorsion) os << " mutorsion " << s.mutorsion;
                if (s.mutorsion!=s.mustorsion) os << " mustorsion " << s.mustorsion;
                if (s.depth) {
                        os << " k2max " << s.k2max 
                                << " kc " << s.kc 
                                << " depth " << s.depth;
                }
                os << " rho " << s.rho 
                        << " dim_particle " << s.dim_particle;
                return os;
        }       
std::istream& operator>> ( std::istream &  is,
CSpecies s 
) [friend]

Reads all species data.

                                                                          {
                is >> s.k >> s.disp >> s.kt >> s.krolling >> s.ktorsion >> s.dispt >> s.disprolling >> s.disptorsion >> s.mu >> s.mus >> s.murolling >> s.musrolling >> s.mutorsion >> s.mustorsion >> s.k2max >> s.kc >> s.depth >> s.rho >> s.dim_particle;
                //s.mus=s.mu;
                return is;
        }

Member Data Documentation

Referenced by CSpecies(), get_rho(), mix(), print(), read(), and set_rho().


The documentation for this class was generated from the following file:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines