|
HG-MD
1
|
This is a particle class. More...
#include <CParticle.h>

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 Particle * | copy () |
| 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 () |
| Particle * | get_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 |
| Particle * | get_PeriodicFromParticle () |
| Vec3D & | get_Position () |
| Mdouble | get_Radius () |
| int | get_Species () |
| int | get_IndSpecies () |
| Vec3D & | get_Velocity () |
| Vec3D & | get_Angle () |
| Vec3D & | get_AngularVelocity () |
| Vec3D & | get_Force () |
| Vec3D & | get_Torque () |
| Vec3D & | get_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) |
| CTangentialSprings & | get_TangentialSprings () |
| CDeltaMaxs & | get_DeltaMaxs () |
Private Attributes | |
| int | HGRID_x |
| Hgrid attributes. | |
| int | HGRID_y |
| int | HGRID_z |
| int | HGRID_Level |
| Cell position in the grid. | |
| Particle * | HGRID_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) | |
| Particle * | periodicFromParticle |
| 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. | |
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
| 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.
{
position = p.position;
velocity = p.velocity;
displacement = p.displacement;
force = p.force;
radius = p.radius;
angle = p.angle;
angularVelocity = p.angularVelocity;
torque = p.torque;
mass = p.get_Mass();
invMass = p.get_InvMass();
inertia = p.get_Inertia();
invInertia = p.get_InvInertia();
HGRID_NextObject = p.HGRID_NextObject;
HGRID_x = p.HGRID_x;
HGRID_y = p.HGRID_y;
HGRID_z = p.HGRID_z;
HGRID_Level = p.HGRID_Level;
tangentialSprings = p.tangentialSprings;
deltaMaxs = p.deltaMaxs;
indSpecies = p.indSpecies;
periodicFromParticle=p.periodicFromParticle;
#ifdef CONSTUCTOR_OUTPUT
cerr << "Particle(Particle &p) finished"<<endl;
#endif
}
| 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
}
| void Particle::accelerate | ( | Vec3D | _new | ) | [inline] |
Referenced by MD::do_integration_after_force_computation(), and MD::do_integration_before_force_computation().
{velocity+=_new;}
| void Particle::add_Force | ( | Vec3D | _new | ) | [inline] |
Referenced by MD::compute_external_forces(), MD::compute_internal_forces(), MD::compute_plastic_internal_forces(), and MD::compute_walls().
{force+=_new;}
| void Particle::add_Torque | ( | Vec3D | _new | ) | [inline] |
Referenced by MD::compute_internal_forces(), MD::compute_plastic_internal_forces(), and MD::compute_walls().
{torque+=_new;}
| void Particle::angularAccelerate | ( | Vec3D | _new | ) | [inline] |
Referenced by MD::do_integration_after_force_computation(), and MD::do_integration_before_force_computation().
{angularVelocity+=_new;}
| 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;}
| Vec3D& Particle::get_AngularVelocity | ( | ) | [inline] |
| CDeltaMaxs& Particle::get_DeltaMaxs | ( | ) | [inline] |
Referenced by MD::compute_plastic_internal_forces().
{return deltaMaxs;}
| Vec3D& Particle::get_Displacement | ( | ) | [inline] |
{return displacement;}
| Vec3D& Particle::get_Force | ( | ) | [inline] |
Referenced by MD::do_integration_after_force_computation(), and MD::do_integration_before_force_computation().
{return force;}
| int Particle::get_HGRID_Level | ( | ) | [inline] |
Referenced by HGRID_2D::CheckCell(), HGRID_3D::CheckCell(), HGRID_2D::CheckObjAgainstGrid(), HGRID_3D::CheckObjAgainstGrid(), HGRID_2D::HGRID_RemoveParticleFromHgrid(), HGRID_3D::HGRID_RemoveParticleFromHgrid(), HGRID_base::HGRID_update_move(), HGRID_2D::HGRID_UpdateParticleInHgrid(), HGRID_3D::HGRID_UpdateParticleInHgrid(), HGRID_2D::TestCell(), and HGRID_3D::TestCell().
{return HGRID_Level;}
| Particle* Particle::get_HGRID_NextObject | ( | ) | [inline] |
Referenced by HGRID_2D::CheckCell(), HGRID_3D::CheckCell(), HGRID_2D::CheckCell_current(), HGRID_3D::CheckCell_current(), HGRID_2D::HGRID_RemoveParticleFromHgrid(), HGRID_3D::HGRID_RemoveParticleFromHgrid(), Chute::IsInsertable(), HGRID_2D::TestCell(), and HGRID_3D::TestCell().
{return HGRID_NextObject;}
| int Particle::get_HGRID_x | ( | ) | [inline] |
| int Particle::get_HGRID_y | ( | ) | [inline] |
| int Particle::get_HGRID_z | ( | ) | [inline] |
Referenced by HGRID_3D::CheckCell(), HGRID_3D::CheckObjAgainstGrid(), HGRID_3D::HGRID_RemoveParticleFromHgrid(), and HGRID_3D::TestCell().
{return HGRID_z;}
| int Particle::get_Index | ( | ) | [inline] |
Referenced by MD::compute_internal_forces(), MD::compute_plastic_internal_forces(), MD::compute_walls(), and HGRID_base::TestObject().
{return index;}
| int Particle::get_IndSpecies | ( | ) | [inline] |
| Mdouble Particle::get_Inertia | ( | ) | const [inline] |
Referenced by Particle().
{return inertia;}
| Mdouble Particle::get_InvInertia | ( | ) | const [inline] |
Referenced by MD::do_integration_after_force_computation(), MD::do_integration_before_force_computation(), and Particle().
{return invInertia;}
| Mdouble Particle::get_InvMass | ( | ) | const [inline] |
Referenced by MD::do_integration_after_force_computation(), MD::do_integration_before_force_computation(), and Particle().
{return invMass;}
| Mdouble Particle::get_KineticEnergy | ( | ) | [inline] |
{return 0.5 * mass * velocity.GetLength2();}
| Mdouble Particle::get_Mass | ( | ) | const [inline] |
Referenced by MD::compute_external_forces(), Chute::get_collision_time(), MD::get_Mass_from_Radius(), MD::get_maximum_velocity(), Chute::get_restitution_coefficient(), Chute::get_SmallestParticle(), Particle(), Chute::set_collision_time_and_restitution_coefficient(), and MD::set_dt().
{return mass;}
| Particle* Particle::get_PeriodicFromParticle | ( | ) | [inline] |
| Vec3D& Particle::get_Position | ( | ) | [inline] |
Referenced by MD::Check_and_Duplicate_Periodic_Particle(), HGRID_2D::CheckObjAgainstGrid(), HGRID_3D::CheckObjAgainstGrid(), Chute::clean_chute(), MD::compute_internal_forces(), MD::compute_plastic_internal_forces(), MD::compute_walls(), Chute::create_bottom(), Chute::create_inflow_particle(), ChuteWithHopper::create_inflow_particle(), CWallPeriodic::distance(), MD::do_integration_after_force_computation(), MD::do_integration_before_force_computation(), CWall::get_distance_and_normal(), CWall::get_distance_and_normal_for_axissymmetric_wall(), HGRID_2D::HGRID_UpdateParticleInHgrid(), HGRID_3D::HGRID_UpdateParticleInHgrid(), ChuteBottom::setup_particles_initial_conditions(), MD::solve(), HGRID_2D::TestObjAgainstGrid(), HGRID_3D::TestObjAgainstGrid(), and HGRID_base::TestObject().
{return position;}
| Mdouble Particle::get_Radius | ( | ) | [inline] |
Referenced by ChuteWithHopper::add_hopper(), MD::Check_and_Duplicate_Periodic_Particle(), MD::compute_internal_forces(), MD::compute_plastic_internal_forces(), MD::compute_walls(), Chute::create_bottom(), Chute::create_inflow_particle(), ChuteWithHopper::create_inflow_particle(), MD::fstat_header(), CWall::get_distance_and_normal(), ParticleHandler::get_LargestParticle(), Chute::get_LargestParticle(), MD::get_maximum_velocity(), Chute::get_radius_of_smallest_particle(), HGrid::InsertParticleToHgrid(), Chute::set_collision_time_and_restitution_coefficient(), ChuteBottom::setup_particles_initial_conditions(), MD::solve(), and HGRID_base::TestObject().
{return radius;}
| int Particle::get_Species | ( | ) | [inline] |
{return indSpecies;}
| CTangentialSprings& Particle::get_TangentialSprings | ( | ) | [inline] |
Referenced by MD::compute_internal_forces(), MD::compute_plastic_internal_forces(), MD::compute_walls(), and Chute::initialize_inflow_particle().
{return tangentialSprings;}
| Vec3D& Particle::get_Torque | ( | ) | [inline] |
Referenced by MD::compute_walls(), MD::do_integration_after_force_computation(), and MD::do_integration_before_force_computation().
{return torque;}
| Vec3D& Particle::get_Velocity | ( | ) | [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] |
Referenced by MD::do_integration_before_force_computation().
{position+=_new;}
| 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().
| 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] |
Referenced by MD::do_integration_after_force_computation(), and MD::do_integration_before_force_computation().
{angle+=_new;}
| void Particle::set_Angle | ( | Vec3D | _new | ) | [inline] |
{angle=_new;}
| void Particle::set_AngularVelocity | ( | Vec3D | _new | ) | [inline] |
{angularVelocity=_new;}
| void Particle::set_Force | ( | Vec3D | _new | ) | [inline] |
{force=_new;}
| void Particle::set_HGRID_Level | ( | int | _new | ) | [inline] |
Referenced by HGrid::InsertParticleToHgrid().
{HGRID_Level=_new;}
| void Particle::set_HGRID_NextObject | ( | Particle * | _new | ) | [inline] |
| void Particle::set_HGRID_x | ( | int | _new | ) | [inline] |
Referenced by HGRID_2D::HGRID_UpdateParticleInHgrid(), and HGRID_3D::HGRID_UpdateParticleInHgrid().
{HGRID_x=_new;}
| void Particle::set_HGRID_y | ( | int | _new | ) | [inline] |
Referenced by HGRID_2D::HGRID_UpdateParticleInHgrid(), and HGRID_3D::HGRID_UpdateParticleInHgrid().
{HGRID_y=_new;}
| void Particle::set_HGRID_z | ( | int | _new | ) | [inline] |
Referenced by HGRID_3D::HGRID_UpdateParticleInHgrid().
{HGRID_z=_new;}
| void Particle::set_Index | ( | int | _new | ) | [inline] |
Referenced by ParticleHandler::copyAndAddParticle().
{index=_new;}
| 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] |
| void Particle::set_periodicFromParticle | ( | Particle * | _new | ) | [inline] |
Referenced by MD::Check_and_Duplicate_Periodic_Particle().
{periodicFromParticle=_new;}
| void Particle::set_Position | ( | Vec3D | _new | ) | [inline] |
| void Particle::set_Radius | ( | Mdouble | _new | ) | [inline] |
Referenced by Chute::create_bottom(), 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(), Chute::set_collision_time_and_restitution_coefficient(), and ChuteBottom::setup_particles_initial_conditions().
{radius=_new;}
| 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);
}
| 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;
}
Vec3D Particle::angle [private] |
Index of the Species of this Particle.
Referenced by Particle().
Vec3D Particle::angularVelocity [private] |
Current angular position.
Referenced by Particle().
CDeltaMaxs Particle::deltaMaxs [private] |
Tangential spring information.
Referenced by Particle().
Vec3D Particle::displacement [private] |
Torque.
Referenced by Particle().
Vec3D Particle::force [private] |
Current particle velocity.
Referenced by Particle().
int Particle::HGRID_Level [private] |
Cell position in the grid.
Referenced by Particle().
Particle* Particle::HGRID_NextObject [private] |
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] |
int Particle::indSpecies [private] |
Pointer to originating Particle.
Referenced by Particle().
Mdouble Particle::inertia [private] |
Inverse Paritcle mass (for computation optimization)
Mdouble Particle::invInertia [private] |
Particle inertia.
Mdouble Particle::invMass [private] |
Particle mass.
Mdouble Particle::mass [private] |
Index of the Particle in the ParticleHandler, used for outputting fstat data.
Particle* Particle::periodicFromParticle [private] |
Particle radius.
Referenced by Particle().
Vec3D Particle::position [private] |
Current angular velocity.
Referenced by Particle().
Mdouble Particle::radius [private] |
Inverse Particle inverse inertia (for computation optimization)
Referenced by Particle().
Displacement (only used in StatisticsVector, StatisticsPoint)
Referenced by Particle().
Vec3D Particle::torque [private] |
Interaction force.
Referenced by Particle().
Vec3D Particle::velocity [private] |
Current particle position.
Referenced by Particle().
1.7.6.1