|
HG-MD
1
|
Chute adds three new effects to the HGrid: the gravity direction can be set using the ChuteAngle variable, a (smooth or rough) bottom wall is created by default, and some basic inflow and outflow routines are added. More...
#include <Chute.h>


Public Member Functions | |
| Chute () | |
| This is the default constructor. All it does is set sensible defaults. | |
| Chute (MD &other) | |
| Copy-constructor for creates an HGRID problem from an existing MD problem. | |
| Chute (HGRID_base &other) | |
| Chute (HGRID_3D &other) | |
| void | constructor () |
| This is the actual constructor; it is called do both constructors above. | |
| void | make_chute_periodic () |
| This makes the chute periodic, in y. | |
| bool | get_IsPeriodic () |
| Get wether the chute is periodic. | |
| void | setup_particles_initial_conditions () |
| initialize particle position, velocity, radius | |
| void | read (std::istream &is) |
| This function reads all chute data. | |
| void | write (std::ostream &os) |
| This function writes all chute data. | |
| void | print (std::ostream &os, bool print_all=false) |
| This function couts all chute data. | |
| void | set_FixedParticleRadius (Mdouble new_) |
| Allows radius of fixed particles to be changed. | |
| Mdouble | get_FixedParticleRadius () |
| Allows radius of fixed particles to be accessed. | |
| void | set_RandomizedBottom (int new_) |
| Changes RandomizedBottom. | |
| int | get_RandomizedBottom () |
| Accesses RandomizedBottom. | |
| void | set_ChuteAngle (Mdouble new_) |
| Sets gravity vector according to chute angle (in degrees) | |
| void | set_ChuteAngle (Mdouble new_, Mdouble gravity) |
| Sets gravity vector according to chute angle (in degrees) | |
| Mdouble | get_ChuteAngle () |
| Gets chute angle (in radians) | |
| Mdouble | get_ChuteAngleDegrees () |
| void | set_max_failed (unsigned int new_) |
| Allows radius of fixed particles to be changed. | |
| unsigned int | get_max_failed () |
| Allows radius of fixed particles to be accessed. | |
| void | set_InflowParticleRadius (Mdouble new_) |
| Allows radius of inflow particles to be changed. | |
| void | set_InflowParticleRadius (Mdouble new_min, Mdouble new_max) |
| Allows radius of inflow particles to be set to a range of values. | |
| void | set_MinInflowParticleRadius (Mdouble new_min) |
| void | set_MaxInflowParticleRadius (Mdouble new_max) |
| Mdouble | get_InflowParticleRadius () |
| Allows radius of inflow particles to be accessed. | |
| Mdouble | get_MinInflowParticleRadius () |
| Allows radius of inflow particles to be accessed. | |
| Mdouble | get_MaxInflowParticleRadius () |
| Allows radius of inflow particles to be accessed. | |
| void | set_InflowHeight (Mdouble new_) |
| Changes inflow height. | |
| Mdouble | get_InflowHeight () |
| Accesses inflow height. | |
| void | set_InflowVelocity (Mdouble new_) |
| Changes inflow velocity. | |
| Mdouble | get_InflowVelocity () |
| Accesses inflow velocity. | |
| void | set_InflowVelocityVariance (Mdouble new_) |
| Changes inflow velocity variance. | |
| Mdouble | get_InflowVelocityVariance () |
| Accesses inflow velocity variance. | |
| void | set_InitialHeight (Mdouble new_) |
| Mdouble | get_InitialHeight () |
| void | set_InitialVelocity (Mdouble new_) |
| Mdouble | get_InitialVelocity () |
| void | set_InitialVelocityVariance (Mdouble new_) |
| Mdouble | get_InitialVelocityVariance () |
| void | set_ChuteWidth (Mdouble new_) |
| Access function that set the width of the chute. | |
| Mdouble | get_ChuteWidth () |
| virtual void | set_ChuteLength (Mdouble new_) |
| Mdouble | get_ChuteLength () |
| int | readNextArgument (unsigned int &i, unsigned int argc, char *argv[]) |
| void | set_collision_time_and_restitution_coefficient (Mdouble tc, Mdouble eps) |
| Sets k, disp such that it matches a given tc and eps for a collision of two inflow particles. | |
| Mdouble | get_collision_time () |
| Calculates collision time of two inflow particles. | |
| Mdouble | get_restitution_coefficient () |
| Calculates restitution coefficient of two inflow particles. | |
| void | set_dt () |
| Sets dt to 1/50-th of the collision time for two particles of mass P. | |
| void | set_dt (Mdouble dt) |
| Sets dt. | |
| Particle * | get_SmallestParticle () |
| Returns the smallest particle (by mass) in the system. | |
| Particle * | get_LargestParticle () |
| Returns the smallest particle (by mass) in the system. | |
| Mdouble | get_radius_of_smallest_particle () |
| Returns the radius of the smallest particle. | |
Protected Member Functions | |
| virtual bool | IsInsertable (Particle &P) |
| here, CheckObjects is called; returns true is the particle should be added | |
| void | add_particle (Particle &P) |
| adds particle to hgrid | |
| void | actions_before_time_step () |
| This is action before the time step is started. | |
| virtual void | add_particles () |
| Here we define the inflow. | |
| void | clean_chute () |
| Here we define the outflow. | |
| void | initialize_inflow_particle () |
| Sets initial values for particles that are created at the inflow. | |
| virtual void | create_inflow_particle () |
| Sets variable values for particles that are created at the inflow. | |
| virtual void | create_bottom () |
| Create the bottom of chute out of particles. | |
| void | cout_time () |
| Couts time. | |
| Mdouble | get_radius_of_largest_particle () |
Protected Attributes | |
| Mdouble | FixedParticleRadius |
| int | RandomizedBottom |
| Mdouble | ChuteAngle |
| Mdouble | MinInflowParticleRadius |
| Mdouble | MaxInflowParticleRadius |
| Mdouble | InflowVelocity |
| Mdouble | InflowVelocityVariance |
| Mdouble | InflowHeight |
| int | max_failed |
| int | num_created |
| Particle | P0 |
Private Attributes | |
| bool | is_periodic |
Chute adds three new effects to the HGrid: the gravity direction can be set using the ChuteAngle variable, a (smooth or rough) bottom wall is created by default, and some basic inflow and outflow routines are added.
| Chute::Chute | ( | ) | [inline] |
This is the default constructor. All it does is set sensible defaults.
References constructor().
{
constructor();
#ifdef CONSTUCTOR_OUTPUT
cerr << "Chute() finished" << endl;
#endif
}
| Chute::Chute | ( | MD & | other | ) | [inline] |
Copy-constructor for creates an HGRID problem from an existing MD problem.
References constructor().
: MD(other), HGRID_3D(other) { constructor(); #ifdef CONSTUCTOR_OUTPUT cerr << "Chute(MD& other) finished" << endl; #endif }
| Chute::Chute | ( | HGRID_base & | other | ) | [inline] |
References constructor().
: MD(other), HGRID_3D(other) { constructor(); #ifdef CONSTUCTOR_OUTPUT cerr << "Chute(HGRID_base& other) finished" << endl; #endif }
| Chute::Chute | ( | HGRID_3D & | other | ) | [inline] |
References constructor().
: MD(other), HGRID_3D(other) { constructor(); #ifdef CONSTUCTOR_OUTPUT cerr << "Chute(HGRID_3D& other) finished" << endl; #endif }
| void Chute::actions_before_time_step | ( | ) | [protected, virtual] |
This is action before the time step is started.
Reimplemented from MD.
Reimplemented in ChuteBottom.
References add_particles(), and clean_chute().
{
add_particles();
clean_chute();
}
| void Chute::add_particle | ( | Particle & | P | ) | [protected] |
adds particle to hgrid
References ParticleHandler::copyAndAddParticle(), MD::get_ParticleHandler(), HGRID_base::grid, HGRID_3D::HGRID_UpdateParticleInHgrid(), and HGrid::InsertParticleToHgrid().
Referenced by IsInsertable().
{
//Puts the particle in the Particle list
get_ParticleHandler().copyAndAddParticle(P);
//This places the particle in this grid
grid->InsertParticleToHgrid(get_ParticleHandler().back());
//This computes where the particle currectly is in the grid
HGRID_UpdateParticleInHgrid(get_ParticleHandler().back());
}
| void Chute::add_particles | ( | ) | [protected, virtual] |
Here we define the inflow.
New particles are created at the inflow, subject to criteria the user can set.
References create_inflow_particle(), IsInsertable(), max_failed, num_created, and P0.
Referenced by actions_before_time_step().
{
int failed = 0;
//try max_failed times to find new insertable particle
while (failed<=max_failed){
create_inflow_particle();
if (IsInsertable(P0)) {
failed = 0;
num_created++;
} else failed++;
};
}
| void Chute::clean_chute | ( | ) | [protected] |
Here we define the outflow.
New particles are destroyed at the outflow, subject to criteria the user can set.
References ParticleHandler::get_NumberOfParticles(), ParticleHandler::get_Particle(), MD::get_ParticleHandler(), Particle::get_Position(), MD::get_xmax(), MD::get_xmin(), MD::removeParticle(), and Vec3D::X.
Referenced by actions_before_time_step().
{
//clean outflow every 100 timesteps
static int count = 0, maxcount = 100;
if (count>maxcount)
{
count = 0;
// delete all outflowing particles
for (unsigned int i=0;i<get_ParticleHandler().get_NumberOfParticles();)
{
if (get_ParticleHandler().get_Particle(i)->get_Position().X>get_xmax()||get_ParticleHandler().get_Particle(i)->get_Position().X<get_xmin())//||get_ParticleHandler().get_Particle(i)->Position.Z+get_ParticleHandler().get_Particle(i)->Radius<zmin)
{
#ifdef DEBUG_OUTPUT_FULL
cout << "erased:" << get_ParticleHandler().get_Particle(i) << endl;
#endif
removeParticle(i);
}
else i++;
}
} else count++;
}
| void Chute::constructor | ( | ) |
This is the actual constructor; it is called do both constructors above.
This is the actually constructor it is called do both constructors above.
Reimplemented from HGRID_3D.
Reimplemented in ChuteWithHopper, and ChuteBottom.
References initialize_inflow_particle(), is_periodic, num_created, set_ChuteAngle(), set_FixedParticleRadius(), set_InflowHeight(), set_InflowParticleRadius(), set_InflowVelocity(), set_InflowVelocityVariance(), set_max_failed(), and set_RandomizedBottom().
Referenced by Chute().
{
is_periodic=false;
set_FixedParticleRadius(0.001);
set_RandomizedBottom(0);
set_ChuteAngle(0.0);
set_max_failed(1);
num_created = 0;
set_InflowParticleRadius(0.001);
set_InflowVelocity(0.1);
set_InflowVelocityVariance(0.0);
set_InflowHeight(0.02);
initialize_inflow_particle();
}
| void Chute::cout_time | ( | ) | [protected, virtual] |
Couts time.
Reimplemented from MD.
References ParticleHandler::get_NumberOfParticles(), MD::get_ParticleHandler(), MD::get_t(), and MD::get_tmax().
{
cout << "t=" << setprecision(3) << left << setw(6) << get_t()
<< ", tmax=" << setprecision(3) << left << setw(6) << get_tmax()
<< ", N=" << setprecision(3) << left << setw(6) << get_ParticleHandler().get_NumberOfParticles()
//~ << ", Nmax=" << setprecision(3) << left << setw(6) << Particles.capacity()
//~ << ", created=" << setprecision(3) << left << setw(6) << num_created
<< endl;
//~ static unsigned int counter=0;
//~ if (++counter>10) {counter=0; cout.flush();}
}
| void Chute::create_bottom | ( | ) | [protected, virtual] |
Create the bottom of chute out of particles.
Creates the bottom of the chute; either smooth, grid-like or random ///.
References ParticleHandler::copyAndAddParticle(), ParticleHandler::end(), get_FixedParticleRadius(), MD::get_ParticleHandler(), Particle::get_Position(), Particle::get_Radius(), RNG::get_RN(), MD::get_xmax(), MD::get_xmin(), MD::get_ymax(), MD::get_ymin(), MD::get_zmin(), HGRID_base::HGRID_actions_before_time_loop(), HGRID_base::HGRID_actions_before_time_step(), IsInsertable(), ChuteBottom::make_rough_bottom(), MD::random, RandomizedBottom, set_InflowParticleRadius(), Particle::set_Position(), Particle::set_Radius(), MD::Walls, Vec3D::X, and Vec3D::Y.
Referenced by setup_particles_initial_conditions(), and ChuteBottom::setup_particles_initial_conditions().
{
if (abs(get_FixedParticleRadius())<1e-12) // smooth bottom
{
//bottom wall
Walls.resize(Walls.size()+1);
Walls.back().set(Vec3D(0.0, 0.0, -1.0), -get_zmin());
}
else //rough bottom
{
// Define standard fixed particle
Particle F0;
F0.set_Radius(get_FixedParticleRadius());
F0.set_Position(Vec3D(0.0,0.0,0.0));
//define grid parameters
Mdouble dx = 2.0 * F0.get_Radius();
Mdouble dy = 2.0 * F0.get_Radius();
int nx = max(1,(int)floor((get_xmax()-get_xmin())/dx));
int ny = max(1,(int)floor((get_ymax()-get_ymin())/dy));
dx = (get_xmax()-get_xmin())/nx; dy = (get_ymax()-get_ymin())/ny;
if (!RandomizedBottom) { // grid-like fixed-particle bottom
for (int i=0; i<nx; i++)
for (int j=0; j<ny; j++)
{
F0.get_Position().X = F0.get_Radius() + dx * i;
F0.get_Position().Y = F0.get_Radius() + dy * j;
get_ParticleHandler().copyAndAddParticle(F0);
}
} else if (RandomizedBottom==1) { // random fixed-particle bottom
cout << "create rough chute bottom, fixed z" << endl;
Mdouble dx = 2.0 * F0.get_Radius();
Mdouble dy = 2.0 * F0.get_Radius();
int nx = max(1,(int)floor((get_xmax()-get_xmin())/dx));
int ny = max(1,(int)floor((get_ymax()-get_ymin())/dy));
dx = (get_xmax()-get_xmin())/nx; dy = (get_ymax()-get_ymin())/ny;
//set_Nmax(2*nx*ny);
//bottom wall
Walls.resize(Walls.size()+1);
Walls.back().set(Vec3D(0.0, 0.0, -1.0), -(get_zmin()-.5*F0.get_Radius()));
//add first particle to initialize HGRID
F0.get_Position().X = random.get_RN(F0.get_Radius(), get_xmax() - F0.get_Radius());
F0.get_Position().Y = random.get_RN(get_ymin()+F0.get_Radius(), get_ymax()-F0.get_Radius());
get_ParticleHandler().copyAndAddParticle(F0);
HGRID_actions_before_time_loop();
HGRID_actions_before_time_step();
//now add more particles
int failed = 0;
while (failed<500)
{
F0.get_Position().X = random.get_RN(F0.get_Radius(), get_xmax()-F0.get_Radius());
F0.get_Position().Y = random.get_RN(get_ymin()+F0.get_Radius(), get_ymax()-F0.get_Radius());
if (IsInsertable(F0)) {
failed = 0;
} else failed++;
}
} else {
//this pointer is the current MD class, so the bottom is create with the particles properties from the MD class
ChuteBottom bottom(*this);
bottom.set_InflowParticleRadius(get_FixedParticleRadius());
bottom.make_rough_bottom(get_ParticleHandler());
cout<<"Starting to destruct ChuteBottom"<<endl;
}
cout<<"Destructed ChuteBottom"<<endl;
//finally, fix particles to the floor
for (vector<Particle*>::iterator it= get_ParticleHandler().begin(); it!=get_ParticleHandler().end(); ++it)
(*it)->fixParticle();
}
}
| void Chute::create_inflow_particle | ( | ) | [protected, virtual] |
Sets variable values for particles that are created at the inflow.
Reimplemented in ChuteWithHopper.
References Particle::compute_particle_mass(), FixedParticleRadius, Particle::get_Position(), Particle::get_Radius(), RNG::get_RN(), Particle::get_Velocity(), MD::get_xmin(), MD::get_ymax(), MD::get_ymin(), MD::get_zmax(), MD::get_zmin(), InflowVelocity, InflowVelocityVariance, MaxInflowParticleRadius, MinInflowParticleRadius, P0, MD::random, Particle::set_Radius(), MD::Species, Vec3D::X, Vec3D::Y, and Vec3D::Z.
Referenced by add_particles().
{
P0.set_Radius(random.get_RN(MinInflowParticleRadius,MaxInflowParticleRadius));
P0.compute_particle_mass(Species);
P0.get_Position().X = get_xmin() + P0.get_Radius();
if ((get_ymax()-get_ymin())==2.0*MaxInflowParticleRadius) P0.get_Position().Y = get_ymin() + P0.get_Radius();
else P0.get_Position().Y = random.get_RN(get_ymin() + P0.get_Radius(), get_ymax() - P0.get_Radius());
P0.get_Position().Z = random.get_RN(get_zmin() + FixedParticleRadius + P0.get_Radius(), get_zmax() - P0.get_Radius());
P0.get_Velocity().X = InflowVelocity * random.get_RN(-InflowVelocityVariance,InflowVelocityVariance) + InflowVelocity;
P0.get_Velocity().Y = InflowVelocity * random.get_RN(-InflowVelocityVariance,InflowVelocityVariance);
P0.get_Velocity().Z = InflowVelocity * random.get_RN(-InflowVelocityVariance,InflowVelocityVariance);
}
| Mdouble Chute::get_ChuteAngle | ( | ) | [inline] |
Gets chute angle (in radians)
References ChuteAngle.
Referenced by ChuteWithHopper::add_hopper(), and ChuteWithHopper::create_inflow_particle().
{return ChuteAngle;}
| Mdouble Chute::get_ChuteAngleDegrees | ( | ) | [inline] |
References ChuteAngle, and constants::pi.
{return ChuteAngle*180./constants::pi;}
| Mdouble Chute::get_ChuteLength | ( | ) | [inline] |
| Mdouble Chute::get_ChuteWidth | ( | ) | [inline] |
References MD::get_ymax().
{return get_ymax();}
Calculates collision time of two inflow particles.
References Particle::compute_particle_mass(), Particle::get_IndSpecies(), Particle::get_Mass(), MaxInflowParticleRadius, MinInflowParticleRadius, P0, Particle::set_Radius(), and MD::Species.
Referenced by ChuteBottom::make_rough_bottom(), and set_dt().
{
P0.set_Radius(0.5*(MinInflowParticleRadius+MaxInflowParticleRadius));
P0.compute_particle_mass(Species);
return MD::get_collision_time(P0.get_Mass(),P0.get_IndSpecies());
}
| Mdouble Chute::get_FixedParticleRadius | ( | ) | [inline] |
Allows radius of fixed particles to be accessed.
References FixedParticleRadius.
Referenced by create_bottom(), and ChuteBottom::setup_particles_initial_conditions().
{return FixedParticleRadius;}
| Mdouble Chute::get_InflowHeight | ( | ) | [inline] |
Accesses inflow height.
References InflowHeight.
Referenced by get_InitialHeight().
{return InflowHeight;}
| Mdouble Chute::get_InflowParticleRadius | ( | ) | [inline] |
Allows radius of inflow particles to be accessed.
References MaxInflowParticleRadius, and MinInflowParticleRadius.
Referenced by ChuteBottom::make_rough_bottom(), and ChuteBottom::setup_particles_initial_conditions().
{return 0.5*(MinInflowParticleRadius+MaxInflowParticleRadius);}
| Mdouble Chute::get_InflowVelocity | ( | ) | [inline] |
Accesses inflow velocity.
References InflowVelocity.
Referenced by get_InitialVelocity().
{return InflowVelocity;}
| Mdouble Chute::get_InflowVelocityVariance | ( | ) | [inline] |
Accesses inflow velocity variance.
References InflowVelocityVariance.
Referenced by get_InitialVelocityVariance().
{return InflowVelocityVariance;}
| Mdouble Chute::get_InitialHeight | ( | ) | [inline] |
References get_InflowHeight().
{
cerr << "WARNING : get_InitialHeight(Mdouble) is a deprecated function, use get_InflowHeight(Mdouble) instead." << endl;
return get_InflowHeight();
}
| Mdouble Chute::get_InitialVelocity | ( | ) | [inline] |
References get_InflowVelocity().
{
cerr << "WARNING : get_InitialVelocity(Mdouble) is a deprecated function, use get_InflowVelocity(Mdouble) instead." << endl;
return get_InflowVelocity();
}
| Mdouble Chute::get_InitialVelocityVariance | ( | ) | [inline] |
References get_InflowVelocityVariance().
{
cerr << "WARNING : get_InitialVelocityVariance(Mdouble) is a deprecated function, use get_InflowVelocityVariance(Mdouble) instead." << endl;
return get_InflowVelocityVariance();
}
| bool Chute::get_IsPeriodic | ( | ) | [inline] |
| Particle * Chute::get_LargestParticle | ( | ) | [virtual] |
Returns the smallest particle (by mass) in the system.
Reimplemented from MD.
References Particle::compute_particle_mass(), ParticleHandler::get_NumberOfParticles(), ParticleHandler::get_Particle(), MD::get_ParticleHandler(), Particle::get_Radius(), MaxInflowParticleRadius, P0, Particle::set_Radius(), and MD::Species.
{
P0.set_Radius(MaxInflowParticleRadius);
P0.compute_particle_mass(Species);
Particle* pP = &P0;
for (unsigned int i=0; i<get_ParticleHandler().get_NumberOfParticles(); i++) {
if (get_ParticleHandler().get_Particle(i)->get_Radius()>pP->get_Radius()) pP = get_ParticleHandler().get_Particle(i);
}
return pP;
}
| unsigned int Chute::get_max_failed | ( | ) | [inline] |
| Mdouble Chute::get_MaxInflowParticleRadius | ( | ) | [inline] |
Allows radius of inflow particles to be accessed.
References MaxInflowParticleRadius.
Referenced by ChuteBottom::setup_particles_initial_conditions().
{return MaxInflowParticleRadius;}
| Mdouble Chute::get_MinInflowParticleRadius | ( | ) | [inline] |
Allows radius of inflow particles to be accessed.
References MinInflowParticleRadius.
Referenced by ChuteBottom::setup_particles_initial_conditions().
{return MinInflowParticleRadius;}
| Mdouble Chute::get_radius_of_largest_particle | ( | ) | [inline, protected] |
References ParticleHandler::end(), FixedParticleRadius, MD::get_ParticleHandler(), and MaxInflowParticleRadius.
{
Mdouble max_rad=max(FixedParticleRadius,MaxInflowParticleRadius);
for (vector<Particle*>::iterator it = get_ParticleHandler().begin(); it!=get_ParticleHandler().end(); ++it)
max_rad=max(max_rad,(*it)->get_Radius());
return max_rad;
}
Returns the radius of the smallest particle.
References FixedParticleRadius, ParticleHandler::get_NumberOfParticles(), ParticleHandler::get_Particle(), MD::get_ParticleHandler(), Particle::get_Radius(), and MinInflowParticleRadius.
{
Mdouble min_rad = MinInflowParticleRadius;
if (FixedParticleRadius) min_rad = min(min_rad,FixedParticleRadius);
for (unsigned int i=0; i<get_ParticleHandler().get_NumberOfParticles(); i++) min_rad=min(min_rad,get_ParticleHandler().get_Particle(i)->get_Radius());
return min_rad;
}
| int Chute::get_RandomizedBottom | ( | ) | [inline] |
Calculates restitution coefficient of two inflow particles.
References Particle::compute_particle_mass(), Particle::get_IndSpecies(), Particle::get_Mass(), MaxInflowParticleRadius, MinInflowParticleRadius, P0, Particle::set_Radius(), and MD::Species.
{
P0.set_Radius(0.5*(MinInflowParticleRadius+MaxInflowParticleRadius));
P0.compute_particle_mass(Species);
return MD::get_restitution_coefficient(P0.get_Mass(),P0.get_IndSpecies());
}
| Particle * Chute::get_SmallestParticle | ( | ) | [virtual] |
Returns the smallest particle (by mass) in the system.
Reimplemented from MD.
References Particle::compute_particle_mass(), Particle::get_Mass(), ParticleHandler::get_NumberOfParticles(), ParticleHandler::get_Particle(), MD::get_ParticleHandler(), MinInflowParticleRadius, P0, Particle::set_Radius(), and MD::Species.
{
P0.set_Radius(MinInflowParticleRadius);
P0.compute_particle_mass(Species);
Particle* pP = &P0;
for (unsigned int i=0; i<get_ParticleHandler().get_NumberOfParticles(); i++) {
if (get_ParticleHandler().get_Particle(i)->get_Mass()<pP->get_Mass()) pP = get_ParticleHandler().get_Particle(i);
}
return pP;
}
| void Chute::initialize_inflow_particle | ( | ) | [protected] |
Sets initial values for particles that are created at the inflow.
References Particle::compute_particle_mass(), Particle::get_Angle(), Particle::get_AngularVelocity(), Particle::get_TangentialSprings(), MinInflowParticleRadius, P0, CTangentialSprings::reset(), Particle::set_Position(), Particle::set_Radius(), Particle::set_Velocity(), Vec3D::set_zero(), and MD::Species.
Referenced by constructor().
{
//Position, maybe also radius and velocity, is reset in create_inflow_particle()
P0.set_Radius(MinInflowParticleRadius);
P0.compute_particle_mass(Species);
P0.get_Angle().set_zero();
P0.get_AngularVelocity().set_zero();
P0.get_TangentialSprings().reset();
P0.set_Position(Vec3D(1e20,1e20,1e20));
P0.set_Velocity(Vec3D(0.0,0.0,0.0));
}
| bool Chute::IsInsertable | ( | Particle & | P | ) | [protected, virtual] |
here, CheckObjects is called; returns true is the particle should be added
todo{Maybe also check if the last particle in a Particular level is removed}
References add_particle(), ParticleHandler::back(), HGrid::ComputeHashBucketIndex(), Particle::get_HGRID_NextObject(), MD::get_ParticleHandler(), HGRID_base::grid, HGrid::objectBucket, ParticleHandler::removeLastParticle(), and HGRID_3D::TestObjAgainstGrid().
Referenced by add_particles(), create_bottom(), and ChuteBottom::setup_particles_initial_conditions().
{
add_particle(P);
if(TestObjAgainstGrid(grid,get_ParticleHandler().back()))
return true;
else
{
Cell cell(get_ParticleHandler().back()->get_HGRID_x(),get_ParticleHandler().back()->get_HGRID_y(), get_ParticleHandler().back()->get_HGRID_z(), get_ParticleHandler().back()->get_HGRID_Level());
int bucket = grid->ComputeHashBucketIndex(cell);
grid->objectBucket[bucket] = get_ParticleHandler().back()->get_HGRID_NextObject();
get_ParticleHandler().removeLastParticle();
return false;
}
}
| void Chute::make_chute_periodic | ( | ) | [inline] |
| void Chute::print | ( | std::ostream & | os, |
| bool | print_all = false |
||
| ) | [virtual] |
This function couts all chute data.
Reimplemented from HGRID_base.
References ChuteAngle, FixedParticleRadius, InflowHeight, InflowVelocity, InflowVelocityVariance, max_failed, MaxInflowParticleRadius, MinInflowParticleRadius, num_created, constants::pi, and RandomizedBottom.
{
HGRID_base::print(os, print_all);
os << " FixedParticleRadius:" << FixedParticleRadius << ", InflowParticleRadius: [" << MinInflowParticleRadius << "," << MaxInflowParticleRadius << "]," << endl
<< " RandomizedBottom:" << RandomizedBottom << ", ChuteAngle:" << ChuteAngle/constants::pi*180. << ", max_failed:" << max_failed << ", num_created:" << num_created << "," << endl
<< " InflowVelocity:" << InflowVelocity << ", InflowVelocityVariance:" << InflowVelocityVariance << ", InflowHeight:" << InflowHeight << endl;
}
| void Chute::read | ( | std::istream & | is | ) | [virtual] |
This function reads all chute data.
Reimplemented from HGRID_base.
Reimplemented in ChuteWithHopper.
References ChuteAngle, FixedParticleRadius, MD::get_restart_version(), InflowHeight, InflowVelocity, InflowVelocityVariance, max_failed, MaxInflowParticleRadius, MinInflowParticleRadius, num_created, constants::pi, and RandomizedBottom.
{
HGRID_base::read(is);
//read out the full line first, so if there is an error it does not affect the read of the next line
string line_string;
getline(is,line_string);
stringstream line (stringstream::in | stringstream::out);
line << line_string;
if (get_restart_version()==1) {
line >> FixedParticleRadius >> RandomizedBottom >> ChuteAngle
>> MinInflowParticleRadius >> MaxInflowParticleRadius >> max_failed >> num_created
>> InflowVelocity >> InflowVelocityVariance >> InflowHeight;
} else {
string dummy;
line >> dummy >> FixedParticleRadius >> dummy >> MinInflowParticleRadius >> dummy >> MaxInflowParticleRadius
>> dummy >> RandomizedBottom >> dummy >> ChuteAngle >> dummy >> max_failed >> dummy >> num_created
>> dummy >> InflowVelocity >> dummy >> InflowVelocityVariance >> dummy >> InflowHeight;
}
//if the Chute Angle is given in degrees, move to radians;
if (ChuteAngle>1.0) ChuteAngle *= constants::pi/180.;
}
| int Chute::readNextArgument | ( | unsigned int & | i, |
| unsigned int | argc, | ||
| char * | argv[] | ||
| ) | [virtual] |
Reimplemented from MD.
Reimplemented in ChuteWithHopper.
References set_ChuteAngle(), set_ChuteLength(), set_ChuteWidth(), set_FixedParticleRadius(), set_InflowHeight(), set_InflowParticleRadius(), MD::set_number_of_saves(), set_RandomizedBottom(), and MD::set_zmax().
{
if (!strcmp(argv[i],"-inflowHeight")) {
set_InflowHeight(atof(argv[i+1]));
set_zmax(atof(argv[i+1]));
} else if (!strcmp(argv[i],"-chuteAngle")) {
set_ChuteAngle(atof(argv[i+1]));
} else if (!strcmp(argv[i],"-chuteLength")) {
set_ChuteLength(atof(argv[i+1]));
} else if (!strcmp(argv[i],"-chuteWidth")) {
set_ChuteWidth(atof(argv[i+1]));
} else if (!strcmp(argv[i],"-number_of_saves")) {
set_number_of_saves(atof(argv[i+1]));
} else if (!strcmp(argv[i],"-fixedParticleRadius")) {
set_FixedParticleRadius(atof(argv[i+1]));
} else if (!strcmp(argv[i],"-inflowParticleRadiusRange")) {
set_InflowParticleRadius(atof(argv[i+1]),atof(argv[i+2]));
i++;
} else if (!strcmp(argv[i],"-inflowParticleRadius")) {
set_InflowParticleRadius(atof(argv[i+1]));
} else if (!strcmp(argv[i],"-randomizedBottom")) {
set_RandomizedBottom(atof(argv[i+1]));
} else return HGRID_3D::readNextArgument(i, argc, argv); //if argv[i] is not found, check the commands in HGRID_3D
return true; //returns true if argv[i] is found
}
| void Chute::set_ChuteAngle | ( | Mdouble | new_ | ) | [inline] |
Sets gravity vector according to chute angle (in degrees)
References MD::get_gravity(), Vec3D::GetLength(), MD::gravity, and set_ChuteAngle().
Referenced by constructor(), ChuteBottom::make_rough_bottom(), readNextArgument(), and set_ChuteAngle().
{Mdouble gravity = get_gravity().GetLength(); if (gravity==0) {cerr<<"WARNING: zero gravity";} set_ChuteAngle(new_, gravity);}
| void Chute::set_ChuteAngle | ( | Mdouble | new_, |
| Mdouble | gravity | ||
| ) | [inline] |
Sets gravity vector according to chute angle (in degrees)
References ChuteAngle, constants::pi, and MD::set_gravity().
{if (new_>=0.0&&new_<=90.0) {ChuteAngle = new_*constants::pi/180.0; set_gravity(Vec3D(sin(ChuteAngle), 0.0, -cos(ChuteAngle))*gravity);} else cerr << "WARNING : Chute angle must be within [0,90]" << endl;}
| virtual void Chute::set_ChuteLength | ( | Mdouble | new_ | ) | [inline, virtual] |
Reimplemented in ChuteWithHopper.
References MD::set_xmax().
Referenced by readNextArgument().
{set_xmax(new_);}
| void Chute::set_ChuteWidth | ( | Mdouble | new_ | ) | [inline] |
Access function that set the width of the chute.
References MD::set_ymax().
Referenced by readNextArgument().
{set_ymax(new_);}
| void Chute::set_collision_time_and_restitution_coefficient | ( | Mdouble | tc, |
| Mdouble | eps | ||
| ) |
Sets k, disp such that it matches a given tc and eps for a collision of two inflow particles.
References Particle::compute_particle_mass(), MD::get_dissipation(), Particle::get_Mass(), Particle::get_Radius(), MD::get_rho(), MaxInflowParticleRadius, MinInflowParticleRadius, P0, constants::pi, MD::set_dissipation(), MD::set_k(), Particle::set_Radius(), MD::Species, and sqr.
Referenced by ChuteBottom::make_rough_bottom().
{
P0.set_Radius(0.5*(MinInflowParticleRadius+MaxInflowParticleRadius));
P0.compute_particle_mass(Species);
set_dissipation(- P0.get_Mass() / tc * log(eps));
set_k(.5 * P0.get_Mass() * (constants::pi*constants::pi/tc/tc + get_dissipation()*get_dissipation()/sqr(P0.get_Mass())));
cout << "collision time and restitution coefficient is set for a particle of radius " << P0.get_Radius() << " and density " << get_rho() << endl;
}
| void Chute::set_dt | ( | ) | [inline] |
Sets dt to 1/50-th of the collision time for two particles of mass P.
Reimplemented from MD.
References get_collision_time(), and set_dt().
Referenced by ChuteBottom::make_rough_bottom(), and set_dt().
{cout<<"Chute set_dt"<<endl;set_dt(.02 * get_collision_time());}
| void Chute::set_dt | ( | Mdouble | dt | ) | [inline] |
| void Chute::set_FixedParticleRadius | ( | Mdouble | new_ | ) | [inline] |
Allows radius of fixed particles to be changed.
References FixedParticleRadius.
Referenced by constructor(), ChuteBottom::make_rough_bottom(), and readNextArgument().
{if (new_ >= 0.0) FixedParticleRadius=new_; else cerr << "WARNING : Fixed particle radius must be greater than or equal to zero" << endl;}
| void Chute::set_InflowHeight | ( | Mdouble | new_ | ) | [inline] |
Changes inflow height.
References InflowHeight, MaxInflowParticleRadius, MinInflowParticleRadius, and MD::set_zmax().
Referenced by constructor(), ChuteBottom::make_rough_bottom(), readNextArgument(), and set_InitialHeight().
{
if (new_ >= MinInflowParticleRadius+MaxInflowParticleRadius) { InflowHeight=new_; set_zmax(1.2*InflowHeight); } else cerr << "WARNING : Inflow height not changed to " << new_ << ", value must be greater than or equal to diameter of inflow particle" << endl;
}
| void Chute::set_InflowParticleRadius | ( | Mdouble | new_ | ) | [inline] |
Allows radius of inflow particles to be changed.
References MaxInflowParticleRadius, and MinInflowParticleRadius.
Referenced by constructor(), create_bottom(), and readNextArgument().
{
if (new_>=0.0) {MinInflowParticleRadius=MaxInflowParticleRadius=new_;} else cerr << "WARNING : Inflow particle must be greater than or equal to zero" << endl;
}
| void Chute::set_InflowParticleRadius | ( | Mdouble | new_min, |
| Mdouble | new_max | ||
| ) | [inline] |
Allows radius of inflow particles to be set to a range of values.
References MaxInflowParticleRadius, and MinInflowParticleRadius.
{
if (new_min>=0.0) {MinInflowParticleRadius=new_min;} else cerr << "WARNING : Min. inflow particle radius must be nonnegative" << endl;
if (new_max>=new_min) {MaxInflowParticleRadius=new_max;} else cerr << "WARNING : Max. inflow particle radius must be >= min. inflow particle radius" << endl;
}
| void Chute::set_InflowVelocity | ( | Mdouble | new_ | ) | [inline] |
Changes inflow velocity.
References InflowVelocity.
Referenced by constructor(), and set_InitialVelocity().
{
if (new_ >= 0.0) InflowVelocity=new_; else cerr << "WARNING : Inflow velocity not changed, value must be greater than or equal to zero" << endl;
}
| void Chute::set_InflowVelocityVariance | ( | Mdouble | new_ | ) | [inline] |
Changes inflow velocity variance.
References InflowVelocityVariance.
Referenced by constructor(), and set_InitialVelocityVariance().
{
if (new_>=0.0&&new_<=1.0) InflowVelocityVariance=new_; else cerr << "WARNING : Inflow velocity variance not changed, value must be within [0,1]" << endl;
}
| void Chute::set_InitialHeight | ( | Mdouble | new_ | ) | [inline] |
References set_InflowHeight().
{
cerr << "WARNING : set_InitialHeight(Mdouble) is a deprecated function, use set_InflowHeight(Mdouble) instead." << endl;
set_InflowHeight(new_);
}
| void Chute::set_InitialVelocity | ( | Mdouble | new_ | ) | [inline] |
References set_InflowVelocity().
{
cerr << "WARNING : set_InitialVelocity(Mdouble) is a deprecated function, use set_InflowVelocity(Mdouble) instead." << endl;
set_InflowVelocity(new_);
}
| void Chute::set_InitialVelocityVariance | ( | Mdouble | new_ | ) | [inline] |
References set_InflowVelocityVariance().
{
cerr << "WARNING : set_InitialVelocityVariance(Mdouble) is a deprecated function, use set_InflowVelocityVariance(Mdouble) instead." << endl;
set_InflowVelocityVariance(new_);
}
| void Chute::set_max_failed | ( | unsigned int | new_ | ) | [inline] |
Allows radius of fixed particles to be changed.
References max_failed.
Referenced by constructor().
{max_failed = new_;}
| void Chute::set_MaxInflowParticleRadius | ( | Mdouble | new_max | ) | [inline] |
References MaxInflowParticleRadius, and MinInflowParticleRadius.
{
if (new_max>=MinInflowParticleRadius) {MaxInflowParticleRadius=new_max;} else cerr << "WARNING : Max. inflow particle radius must be >= min. inflow particle radius" << endl;
}
| void Chute::set_MinInflowParticleRadius | ( | Mdouble | new_min | ) | [inline] |
References MaxInflowParticleRadius, and MinInflowParticleRadius.
{
if (new_min<=MaxInflowParticleRadius) {MinInflowParticleRadius=new_min;} else cerr << "WARNING : Max. inflow particle radius must be >= min. inflow particle radius" << endl;
}
| void Chute::set_RandomizedBottom | ( | int | new_ | ) | [inline] |
Changes RandomizedBottom.
References RandomizedBottom.
Referenced by constructor(), ChuteBottom::make_rough_bottom(), and readNextArgument().
{RandomizedBottom = new_;}
| void Chute::setup_particles_initial_conditions | ( | ) | [virtual] |
initialize particle position, velocity, radius
This initially set up the particles///.
Reimplemented from MD.
Reimplemented in ChuteBottom, and ChuteWithHopper.
References create_bottom(), MD::get_ymax(), MD::get_ymin(), is_periodic, MD::set_NWall(), MD::set_NWallPeriodic(), MD::Walls, and MD::WallsPeriodic.
{
//set side walls - solid if not a periodic
if (is_periodic)
{
set_NWallPeriodic(1);
set_NWall(0);
WallsPeriodic[0].set(Vec3D( 0.0, 1.0, 0.0), get_ymin(), get_ymax());
}
else
{
set_NWallPeriodic(0);
set_NWall(2);
Walls[0].set(Vec3D( 0.0,-1.0, 0.0), -get_ymin());
Walls[1].set(Vec3D( 0.0, 1.0, 0.0), get_ymax());
}
//creates the bottom of the chute
create_bottom();
}
| void Chute::write | ( | std::ostream & | os | ) | [virtual] |
This function writes all chute data.
Reimplemented from HGRID_base.
Reimplemented in ChuteWithHopper.
References ChuteAngle, FixedParticleRadius, InflowHeight, InflowVelocity, InflowVelocityVariance, max_failed, MaxInflowParticleRadius, MinInflowParticleRadius, num_created, constants::pi, and RandomizedBottom.
{
HGRID_base::write(os);
os<< "FixedParticleRadius " << FixedParticleRadius
<< " MinInflowParticleRadius " << MinInflowParticleRadius
<< " MaxInflowParticleRadius " << MaxInflowParticleRadius
<< endl
<< "RandomizedBottom " << RandomizedBottom
<< " ChuteAngle " << ChuteAngle/constants::pi*180.
<< " max_failed " << max_failed
<< " num_created " << num_created
<< endl
<< "InflowVelocity " << InflowVelocity
<< " InflowVelocityVariance " << InflowVelocityVariance
<< " InflowHeight " << InflowHeight
<< endl;
}
Mdouble Chute::ChuteAngle [protected] |
Mdouble Chute::FixedParticleRadius [protected] |
Mdouble Chute::InflowHeight [protected] |
Referenced by get_InflowHeight(), print(), read(), set_InflowHeight(), and write().
Mdouble Chute::InflowVelocity [protected] |
Referenced by create_inflow_particle(), get_InflowVelocity(), print(), read(), set_InflowVelocity(), and write().
Mdouble Chute::InflowVelocityVariance [protected] |
Referenced by create_inflow_particle(), get_InflowVelocityVariance(), print(), read(), set_InflowVelocityVariance(), and write().
bool Chute::is_periodic [private] |
Referenced by constructor(), get_IsPeriodic(), make_chute_periodic(), and setup_particles_initial_conditions().
int Chute::max_failed [protected] |
Referenced by add_particles(), get_max_failed(), print(), read(), set_max_failed(), ChuteBottom::setup_particles_initial_conditions(), and write().
Mdouble Chute::MaxInflowParticleRadius [protected] |
Referenced by create_inflow_particle(), ChuteWithHopper::create_inflow_particle(), get_collision_time(), get_InflowParticleRadius(), get_LargestParticle(), get_MaxInflowParticleRadius(), get_radius_of_largest_particle(), get_restitution_coefficient(), ChuteBottom::make_rough_bottom(), print(), read(), set_collision_time_and_restitution_coefficient(), set_InflowHeight(), set_InflowParticleRadius(), set_MaxInflowParticleRadius(), set_MinInflowParticleRadius(), and write().
Mdouble Chute::MinInflowParticleRadius [protected] |
Referenced by create_inflow_particle(), ChuteWithHopper::create_inflow_particle(), get_collision_time(), get_InflowParticleRadius(), get_MinInflowParticleRadius(), get_radius_of_smallest_particle(), get_restitution_coefficient(), get_SmallestParticle(), initialize_inflow_particle(), print(), read(), set_collision_time_and_restitution_coefficient(), set_InflowHeight(), set_InflowParticleRadius(), set_MaxInflowParticleRadius(), set_MinInflowParticleRadius(), and write().
int Chute::num_created [protected] |
Referenced by add_particles(), constructor(), print(), read(), ChuteBottom::setup_particles_initial_conditions(), and write().
Referenced by ChuteWithHopper::add_hopper(), add_particles(), create_inflow_particle(), ChuteWithHopper::create_inflow_particle(), get_collision_time(), get_LargestParticle(), get_restitution_coefficient(), get_SmallestParticle(), initialize_inflow_particle(), set_collision_time_and_restitution_coefficient(), and ChuteBottom::setup_particles_initial_conditions().
int Chute::RandomizedBottom [protected] |
Referenced by create_bottom(), get_RandomizedBottom(), print(), read(), set_RandomizedBottom(), and write().
1.7.6.1