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

This is a particle class. More...

#include <CParticle.h>

Collaboration diagram for Particle:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 Particle ()
 Basic Particle contructors, creates an Particle at (0,0,0) with radius, mass and inertia equal to 1.
 Particle (const Particle &p)
 Particle copy contstructor, which accepts as input a reference to a Particle. It creates a copy of this Particle and all it's information. Usually it is better to use the copy() function for polymorfism.
virtual ~Particle ()
 Particle destructor, needs to be implemented and checked if it removes tangential spring information.
virtual Particlecopy ()
 Particle copy method. It calls to copy contrustor of this Particle, usefull for polymorfism.
Mdouble get_Volume (vector< CSpecies > &Species) const
 Get Particle volume function, which required a reference to the Species vector. It returns the volume of the Particle.
void fixParticle ()
 Fix Particle function. It fixes a Particle by setting its inverse mass and inertia and velocities to zero.
bool isFixed ()
 Is fixed Particle function. It returns wether a Particle is fixed or not, by cheking its inverse Mass.
void unfix (vector< CSpecies > &Species)
 Unfix Particle function, which required a reference to the Species vector. It un fixes a Particle by compyting the Particles mass and inertia.
void compute_particle_mass (vector< CSpecies > &Species)
 Compute Particle mass function, which required a reference to the Species vector. It copmuters the Particles mass, Inertia and the inverses.
void print (std::ostream &os)
 Particle print function, which accepts an os stringstream as input. It prints human readable Particle information to the stringstream.
void print_HGRID (std::ostream &os)
int get_HGRID_Level ()
Particleget_HGRID_NextObject ()
int get_HGRID_x ()
int get_HGRID_y ()
int get_HGRID_z ()
int get_Index ()
Mdouble get_Inertia () const
Mdouble get_InvInertia () const
Mdouble get_InvMass () const
Mdouble get_KineticEnergy ()
Mdouble get_Mass () const
Particleget_PeriodicFromParticle ()
Vec3Dget_Position ()
Mdouble get_Radius ()
int get_Species ()
int get_IndSpecies ()
Vec3Dget_Velocity ()
Vec3Dget_Angle ()
Vec3Dget_AngularVelocity ()
Vec3Dget_Force ()
Vec3Dget_Torque ()
Vec3Dget_Displacement ()
void set_inertia (Mdouble _new)
void set_periodicFromParticle (Particle *_new)
void set_species (Mdouble _new)
void set_Index (int _new)
void set_HGRID_x (int _new)
void set_HGRID_y (int _new)
void set_HGRID_z (int _new)
void set_HGRID_Level (int _new)
void set_HGRID_NextObject (Particle *_new)
void set_Radius (Mdouble _new)
void set_IndSpecies (int _new)
void set_Mass (Mdouble _new)
void set_Velocity (Vec3D _new)
void set_Position (Vec3D _new)
void set_Angle (Vec3D _new)
void set_AngularVelocity (Vec3D _new)
void set_Force (Vec3D _new)
void set_Torque (Vec3D _new)
void move (Vec3D _new)
void accelerate (Vec3D _new)
void rotate (Vec3D _new)
void angularAccelerate (Vec3D _new)
void add_Force (Vec3D _new)
void add_Torque (Vec3D _new)
CTangentialSpringsget_TangentialSprings ()
CDeltaMaxsget_DeltaMaxs ()

Private Attributes

int HGRID_x
 Hgrid attributes.
int HGRID_y
int HGRID_z
int HGRID_Level
 Cell position in the grid.
ParticleHGRID_NextObject
 Grid level for the object.
int index
 Pointer to the next Particle in the same HGrid cell.
Mdouble mass
 Index of the Particle in the ParticleHandler, used for outputting fstat data.
Mdouble invMass
 Particle mass.
Mdouble inertia
 Inverse Paritcle mass (for computation optimization)
Mdouble invInertia
 Particle inertia.
Mdouble radius
 Inverse Particle inverse inertia (for computation optimization)
ParticleperiodicFromParticle
 Particle radius.
int indSpecies
 Pointer to originating Particle.
Vec3D angle
 Index of the Species of this Particle.
Vec3D angularVelocity
 Current angular position.
Vec3D position
 Current angular velocity.
Vec3D velocity
 Current particle position.
Vec3D force
 Current particle velocity.
Vec3D torque
 Interaction force.
Vec3D displacement
 Torque.
CTangentialSprings tangentialSprings
 Displacement (only used in StatisticsVector, StatisticsPoint)
CDeltaMaxs deltaMaxs
 Tangential spring information.

Friends

std::ostream & operator<< (std::ostream &os, const Particle &p)
 Particle << opperator, writes Particle information to the stringstream.
std::istream & operator>> (std::istream &is, Particle &p)
 Particle >> opperator, reads Particle information from the stringstream.

Detailed Description

This is a particle class.

This code is orginally from Vitaliy. It also comes with a required Vector class Modifications, 20:8/9:2009 Included the Vector class here rather that at the level above and hence, added inclusion guards


Constructor & Destructor Documentation

Particle::Particle ( ) [inline]

Basic Particle contructors, creates an Particle at (0,0,0) with radius, mass and inertia equal to 1.

                   {
                position.set_zero();
                velocity.set_zero();
                displacement.set_zero();
                force.set_zero();
                radius=1.0;
                angle.set_zero();
                angularVelocity.set_zero();
                torque.set_zero();
                mass = invMass = 1.0;
                inertia = invInertia = 1.0;
                HGRID_NextObject=NULL;
                
                tangentialSprings.reset();
                deltaMaxs.reset();
                indSpecies = 0;
                periodicFromParticle=NULL;
                
                #ifdef CONSTUCTOR_OUTPUT
                        cerr << "Particle() finished"<<endl;
                #endif          
        }
Particle::Particle ( const Particle p) [inline]

Particle copy contstructor, which accepts as input a reference to a Particle. It creates a copy of this Particle and all it's information. Usually it is better to use the copy() function for polymorfism.

References angle, angularVelocity, deltaMaxs, displacement, force, get_Inertia(), get_InvInertia(), get_InvMass(), get_Mass(), HGRID_Level, HGRID_NextObject, HGRID_x, HGRID_y, HGRID_z, indSpecies, periodicFromParticle, position, radius, tangentialSprings, torque, and velocity.

virtual Particle::~Particle ( ) [inline, virtual]

Particle destructor, needs to be implemented and checked if it removes tangential spring information.

        {
                #ifdef CONSTUCTOR_OUTPUT
                        cerr << "virtual ~Particle finished"<<endl;
                #endif  
        }

Member Function Documentation

void Particle::accelerate ( Vec3D  _new) [inline]
void Particle::add_Force ( Vec3D  _new) [inline]
void Particle::add_Torque ( Vec3D  _new) [inline]
void Particle::angularAccelerate ( Vec3D  _new) [inline]
void Particle::compute_particle_mass ( vector< CSpecies > &  Species) [inline]

Compute Particle mass function, which required a reference to the Species vector. It copmuters the Particles mass, Inertia and the inverses.

References constants::pi, and sqr.

Referenced by Chute::create_inflow_particle(), ChuteWithHopper::create_inflow_particle(), Chute::get_collision_time(), Chute::get_LargestParticle(), MD::get_Mass_from_Radius(), Chute::get_restitution_coefficient(), Chute::get_SmallestParticle(), Chute::initialize_inflow_particle(), MD::read(), MD::read_v1(), MD::read_v2(), Chute::set_collision_time_and_restitution_coefficient(), and ChuteBottom::setup_particles_initial_conditions().

                                                              {
                if (!isFixed()) {
                        switch(Species[indSpecies].dim_particle)
                        {
                                case 3 :
                                { 
                                        set_Mass(4.0/3.0 * constants::pi * radius * radius * radius * Species[indSpecies].rho); 
                                        set_inertia(.4 * get_Mass() * sqr(radius));
                                        break;
                                }
                                case 2 :
                                { 
                                        set_Mass(constants::pi * radius * radius * Species[indSpecies].rho);                              
                                        set_inertia(.5 * get_Mass() * sqr(radius));
                                        break;
                                }
                                case 1 :
                                {
                                        set_Mass(2.0 * radius * Species[indSpecies].rho);                                     
                                        set_inertia(0.0);
                                        break;
                                }
                                default :
                                {
                                        cerr<<"In compute_particle_mass(vector<CSpecies>& Species) the dimension of the particle is not set"<<endl;
                                        exit(-1);
                                }
                        }
                }
        }
virtual Particle* Particle::copy ( ) [inline, virtual]

Particle copy method. It calls to copy contrustor of this Particle, usefull for polymorfism.

Referenced by ParticleHandler::copyAndAddParticle().

        {
                return new Particle(*this);
        }
void Particle::fixParticle ( ) [inline]

Fix Particle function. It fixes a Particle by setting its inverse mass and inertia and velocities to zero.

Referenced by MD::read_next_from_data_file().

                          {
                mass=1e20; invMass=0.0; velocity.set_zero(); inertia=1e20; invInertia=0.0; angularVelocity.set_zero();
#ifdef USE_SIMPLE_VERLET_INTEGRATION
                position = PrevPosition;
#endif  
        }
Vec3D& Particle::get_Angle ( ) [inline]

Referenced by Chute::initialize_inflow_particle().

{return angle;}
{return displacement;}
Vec3D& Particle::get_Force ( ) [inline]
int Particle::get_HGRID_Level ( ) [inline]
int Particle::get_HGRID_x ( ) [inline]
int Particle::get_HGRID_y ( ) [inline]
int Particle::get_HGRID_z ( ) [inline]
int Particle::get_Index ( ) [inline]
int Particle::get_IndSpecies ( ) [inline]
Mdouble Particle::get_Inertia ( ) const [inline]

Referenced by Particle().

{return inertia;}
Mdouble Particle::get_InvInertia ( ) const [inline]
Mdouble Particle::get_InvMass ( ) const [inline]
{return 0.5 * mass * velocity.GetLength2();}
Mdouble Particle::get_Mass ( ) const [inline]
int Particle::get_Species ( ) [inline]
{return indSpecies;}
Vec3D& Particle::get_Torque ( ) [inline]
Mdouble Particle::get_Volume ( vector< CSpecies > &  Species) const [inline]

Get Particle volume function, which required a reference to the Species vector. It returns the volume of the Particle.

References constants::pi.

        {
                switch(Species[indSpecies].dim_particle)
                {
                        case 3 :
                                return (4.0/3.0 * constants::pi * radius * radius * radius); 
                        case 2 :
                                return (constants::pi * radius * radius);                              
                        case 1 :
                                return (2.0 * radius);                                     
                        default :
                        {
                                cerr<<"In get_Volume(vector<CSpecies>& Species) the dimension of the particle is not set"<<endl;
                                exit(-1);
                        }
                }
                
        } 
bool Particle::isFixed ( ) [inline]

Is fixed Particle function. It returns wether a Particle is fixed or not, by cheking its inverse Mass.

Referenced by MD::compute_external_forces(), MD::compute_internal_forces(), MD::compute_plastic_internal_forces(), and MD::compute_walls().

{return (invMass==0.0);}
void Particle::move ( Vec3D  _new) [inline]
void Particle::print ( std::ostream &  os) [inline]

Particle print function, which accepts an os stringstream as input. It prints human readable Particle information to the stringstream.

Referenced by MD::print().

                                   {
                os << "Particle( Position:" << position 
                        << ", Velocity:" << velocity 
                        << ", Radius:" << radius 
                        << ", Mass:" << mass 
         //<< ", InvMass:" << invMass
         //<< ", Angle:" << Angle 
                        << ", AngularVelocity:" << angularVelocity 
                        << ", Inertia:" << inertia 
                        << ")";
        }
void Particle::print_HGRID ( std::ostream &  os) [inline]
                                         {
                os << "Particle( HGRID_Level:" << HGRID_Level 
                        << ", HGRID_x:" << HGRID_x 
                        << ", HGRID_y:" << HGRID_y 
                        << ", HGRID_z:" << HGRID_z 
                        << ")";
        }       
void Particle::rotate ( Vec3D  _new) [inline]
void Particle::set_Angle ( Vec3D  _new) [inline]
{angle=_new;}
void Particle::set_AngularVelocity ( Vec3D  _new) [inline]
void Particle::set_Force ( Vec3D  _new) [inline]
{force=_new;}
void Particle::set_HGRID_Level ( int  _new) [inline]
void Particle::set_HGRID_NextObject ( Particle _new) [inline]
void Particle::set_HGRID_x ( int  _new) [inline]
void Particle::set_HGRID_y ( int  _new) [inline]
void Particle::set_HGRID_z ( int  _new) [inline]
void Particle::set_Index ( int  _new) [inline]
void Particle::set_IndSpecies ( int  _new) [inline]

Referenced by MD::get_Mass_from_Radius().

{indSpecies=_new;}
void Particle::set_inertia ( Mdouble  _new) [inline]
{if (_new>=0){inertia=_new; invInertia=1.0/_new;} else { cerr << "Error in set_inertia ("<<_new<<")"<< endl; exit(-1); }}
void Particle::set_Mass ( Mdouble  _new) [inline]
{if(_new>=0.0) {if(invMass){mass=_new;invMass=1.0/_new;}} else { cerr << "Error in set_Mass(" << _new << ")" << endl; exit(-1); }} //InvMass=0 is a flag for a fixed particle
void Particle::set_periodicFromParticle ( Particle _new) [inline]
void Particle::set_Position ( Vec3D  _new) [inline]
void Particle::set_Radius ( Mdouble  _new) [inline]
void Particle::set_species ( Mdouble  _new) [inline]
{indSpecies=_new;}
void Particle::set_Torque ( Vec3D  _new) [inline]

Referenced by MD::compute_walls().

{torque=_new;}
void Particle::set_Velocity ( Vec3D  _new) [inline]
void Particle::unfix ( vector< CSpecies > &  Species) [inline]

Unfix Particle function, which required a reference to the Species vector. It un fixes a Particle by compyting the Particles mass and inertia.

                                             {
                invMass=1.0;
                compute_particle_mass(Species);
        }

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  os,
const Particle p 
) [friend]

Particle << opperator, writes Particle information to the stringstream.

todo: add plastic history parameter

        {
                os << p.position << " " << p.velocity << " " << p.radius << " " << p.angle << " " << p.angularVelocity << " " << p.invMass << " " << p.invInertia << " " << p.tangentialSprings;
                return os;
        }       
std::istream& operator>> ( std::istream &  is,
Particle p 
) [friend]

Particle >> opperator, reads Particle information from the stringstream.

        {
                is >> p.position >> p.velocity >> p.radius >> p.angle >> p.angularVelocity >> p.invMass >> p.invInertia >> p.tangentialSprings; //add this for the new restart
                if (p.invMass) p.mass = 1.0/p.invMass; else p.mass = 1e20; 
                if (p.invInertia) p.inertia = 1.0/p.invInertia; else p.inertia = 1e20; 
                return is;
        }

Member Data Documentation

Index of the Species of this Particle.

Referenced by Particle().

Current angular position.

Referenced by Particle().

Tangential spring information.

Referenced by Particle().

Torque.

Referenced by Particle().

Current particle velocity.

Referenced by Particle().

int Particle::HGRID_Level [private]

Cell position in the grid.

Referenced by Particle().

Grid level for the object.

Referenced by Particle().

int Particle::HGRID_x [private]

Hgrid attributes.

Referenced by Particle().

int Particle::HGRID_y [private]

Referenced by Particle().

int Particle::HGRID_z [private]

Referenced by Particle().

int Particle::index [private]

Pointer to the next Particle in the same HGrid cell.

Particle attributes

int Particle::indSpecies [private]

Pointer to originating Particle.

Referenced by Particle().

Inverse Paritcle mass (for computation optimization)

Particle inertia.

Particle mass.

Index of the Particle in the ParticleHandler, used for outputting fstat data.

Particle radius.

Referenced by Particle().

Current angular velocity.

Referenced by Particle().

Inverse Particle inverse inertia (for computation optimization)

Referenced by Particle().

Displacement (only used in StatisticsVector, StatisticsPoint)

Referenced by Particle().

Interaction force.

Referenced by Particle().

Current particle position.

Referenced by Particle().


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