HG-MD  1
Public Member Functions | Protected Attributes | Private Attributes
ChuteWithHopper Class Reference

ChuteWithHopper has a hopper as inflow. More...

#include <ChuteWithHopper.h>

Inheritance diagram for ChuteWithHopper:
Inheritance graph
[legend]
Collaboration diagram for ChuteWithHopper:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 ChuteWithHopper (Chute &other)
 This is a copy constructor for Chute problems.
 ChuteWithHopper (HGRID_3D &other)
 ChuteWithHopper (HGRID_base &other)
 ChuteWithHopper (MD &other)
 ChuteWithHopper ()
 This is the default constructor.
void constructor ()
 This is the actually constructor, get called by all constructors above.
void set_HopperFillPercentage (Mdouble new_fill)
virtual void setup_particles_initial_conditions ()
 initialize particle position, velocity, radius
void add_hopper ()
 This create the hopper on top of the chute, see diagram in class description for details of the points.
virtual void create_inflow_particle ()
 This creates an inflow particle in the top 50% of the hopper i.e.
void set_Hopper (Mdouble ExitLength, Mdouble ExitHeight, Mdouble Angle, Mdouble Length)
Mdouble get_MaximumVelocityInducedByGravity ()
 Allows chute length to be accessed.
Mdouble get_ChuteLength ()
 Allows chute length to be accessed.
void set_ChuteLength (Mdouble new_)
 Allows chute length to be changed.
void set_centerHopper (bool new_)
void set_lowerFillHeight (Mdouble new_)
void set_shift (Mdouble new_)
virtual void read (std::istream &is)
 This function reads all chute data.
virtual void write (std::ostream &os)
 This function writes all chute data.
virtual void print (std::ostream &os)
void lift_hopper (Mdouble distance)
 This lifts the hopper above the plane of the chute.
Mdouble get_lift_hopper ()
void set_hopper_dim (Mdouble new_hopper_dim)
void set_align_base (bool new_align)
int readNextArgument (unsigned int &i, unsigned int argc, char *argv[])

Protected Attributes

Mdouble HopperLength
 Dimension of the hopper in vertical direction.
Mdouble HopperHeight
 Dimension of the hopper in horizontal direction.
Mdouble HopperAngle
 Angle between the two pieces of the hopper walls.
Mdouble HopperExitLength
 Dimension of the hopper exit in vertical direction.
Mdouble HopperExitHeight
 Dimension of the hopper exit in vertical direction.
Mdouble shift
 The x position where the Chute starts (defined as the beginning of the hopper)
Mdouble lowerFillHeight
 Relative height (in [0,1)) above which teh hopper is replenished with new particles.
bool centerHopper
 If theis flag is set, the hopper will be constructed in the xy-center of the domain, and not next to the xmin-domain boundary; by default off.

Private Attributes

Mdouble lift
 This is the vertical distance the chute is lifted above the plane.
unsigned int hopper_dim
 This is the dimension of the hopper, my default it is one dimensional and hence does not have side wall.
bool align_base
 This is the flag, which sets if the chute bottom is aligned with the hopper, by default it is.
Mdouble fill_percent
 This is which percentage of the hopper is used for creating new partices;.

Detailed Description

ChuteWithHopper has a hopper as inflow.

The hopper has two parts as follows to create the finite hopper walls, we take vector between two wall points in xz-plane, then rotate clockwise and make unit length.

hopper.jpg
Sketch of the hopper

A,B,C denote three points on the left and right hopper walls which are used to construct the hopper shift denotes the space by which the chute has to be shifted to the right such that the hopper is in the domain Note the wall direction has to be set seperately either period of walls.


Constructor & Destructor Documentation

ChuteWithHopper::ChuteWithHopper ( Chute other) [inline]

This is a copy constructor for Chute problems.

Bug:
This copy construct is untested

References constructor().

:               MD(other), Chute(other) {constructor();}

References constructor().

:               MD(other), Chute(other) {constructor();}

References constructor().

:       MD(other), Chute(other) {constructor();}
ChuteWithHopper::ChuteWithHopper ( MD other) [inline]

References constructor().

:                       MD(other), Chute(other) {constructor();}

This is the default constructor.

References constructor().


Member Function Documentation

void ChuteWithHopper::add_hopper ( ) [inline]

This create the hopper on top of the chute, see diagram in class description for details of the points.

References centerHopper, Chute::ChuteAngle, ParticleHandler::end(), Chute::get_ChuteAngle(), MD::get_NWall(), MD::get_ParticleHandler(), Particle::get_Radius(), MD::get_ymax(), MD::get_ymin(), Vec3D::GetLength2, hopper_dim, HopperAngle, HopperExitHeight, HopperExitLength, HopperHeight, HopperLength, lift, Chute::P0, MD::set_NWall(), set_shift(), MD::set_zmax(), shift, MD::Walls, Vec3D::X, and Vec3D::Z.

Referenced by setup_particles_initial_conditions().

                          {
                //hopper walls
                int n = get_NWall();
                set_NWall(n+2);
                
                //to create the finite hopper walls, we take vector between two wall points in xz-plane, then rotate clockwise and make unit length
                // A\       /A  
                //   \     /       A,B,C denote three points on the left and right hopper walls which are used to construct the hopper
                //    \   /        shift denotes the space by which the chute has to be shifted to the right such that the hopper is in the domain
                //   B|   |B
                //    |   |
                //    |   |
                //   C|   |C
                
                Vec3D A, B, C, temp, normal;
                Mdouble s = sin(get_ChuteAngle());
                Mdouble c = cos(get_ChuteAngle());
                Mdouble HopperLowestPoint = HopperExitHeight - HopperExitLength * tan(ChuteAngle);
                HopperHeight = HopperLowestPoint + 1.1 * 0.5*(HopperLength+HopperExitLength) / tan(HopperAngle);
                Mdouble HopperCornerHeight = HopperHeight - 0.5*(HopperLength-HopperExitLength) / tan(HopperAngle);
                if (HopperCornerHeight<=0.0) { HopperHeight += -HopperCornerHeight + P0.get_Radius(); HopperCornerHeight = P0.get_Radius(); }
                
                //first we create the left hopper wall
                
                //coordinates of A,B,C in (vertical parallel to flow,vertical normal to flow, horizontal) direction
                A = Vec3D(-0.5*(HopperLength-HopperExitLength), 0.0, HopperHeight);
                B = Vec3D(0.0, 0.0, HopperCornerHeight);
                C = Vec3D(0.0, 0.0, 0.0);
                
                
                
                //now rotate the coordinates of A,B,C to be in (x,y,z) direction                
                A = Vec3D(c*A.X-s*A.Z, 0.0, s*A.X+c*A.Z);
                B = Vec3D(c*B.X-s*B.Z, 0.0, s*B.X+c*B.Z);
                C = Vec3D(c*C.X-s*C.Z, 0.0, s*C.X+c*C.Z);
                // the position of A determines shift and zmax
                if (centerHopper) set_shift(-A.X+40);
                else set_shift(-A.X);
                set_zmax(A.Z);
                A.X +=shift;
                B.X +=shift;
                C.X +=shift;
                
                //This lifts the hopper a distance above the chute
                A.Z+=lift;
                B.Z+=lift;
                C.Z+=lift;
                
                //create a finite wall from B to A and from C to B
                temp = B-A;
                normal  = Vec3D(temp.Z,0.0,-temp.X) / sqrt(temp.GetLength2());
                Walls[n].add_finite_wall(normal, Dot(normal,A));
                temp = C-B;
                normal  = Vec3D(temp.Z,0.0,-temp.X) / sqrt(temp.GetLength2());
                Walls[n].add_finite_wall(normal, Dot(normal,B));
                temp = A-C;
                normal = Vec3D(temp.Z,0.0,-temp.X)/sqrt(temp.GetLength2());
                Walls[n].add_finite_wall(normal,Dot(normal,C));
                
                
                
                //next, do the same for the right wall
                A = Vec3D(HopperLength-0.5*(HopperLength-HopperExitLength), 0.0, HopperHeight);
                B = Vec3D(0.5*(HopperLength+HopperExitLength)-0.5*(HopperLength-HopperExitLength), 0.0, HopperCornerHeight);
                C = Vec3D(0.5*(HopperLength+HopperExitLength)-0.5*(HopperLength-HopperExitLength), 0.0, HopperLowestPoint);
                
                
                
                
                //This rotates the right points
                A = Vec3D(c*A.X-s*A.Z+shift, 0.0, s*A.X+c*A.Z);
                B = Vec3D(c*B.X-s*B.Z+shift, 0.0, s*B.X+c*B.Z);
                C = Vec3D(c*C.X-s*C.Z+shift, 0.0, s*C.X+c*C.Z);
                
                //This lifts the hopper a distance above the chute
                A.Z+=lift;
                B.Z+=lift;
                C.Z+=lift;
                
                temp = A-B;
                normal  = Vec3D(temp.Z,0.0,-temp.X) / sqrt(temp.GetLength2());
                Walls[n+1].add_finite_wall(normal, Dot(normal,A));
                temp = B-C;
                normal  = Vec3D(temp.Z,0.0,-temp.X) / sqrt(temp.GetLength2());
                Walls[n+1].add_finite_wall(normal, Dot(normal,B));
                temp = C-A;
                normal  = Vec3D(temp.Z,0.0,-temp.X) / sqrt(temp.GetLength2());
                Walls[n+1].add_finite_wall(normal, Dot(normal,C));
                
                set_zmax(A.Z);
                
                if (hopper_dim == 2)
                        {
                                
                                set_NWall(n+4);
                                
                                
                                //coordinates of A,B,C in (vertical parallel to flow,vertical normal to flow, horizontal) direction
                                A = Vec3D(0.0, (get_ymax()-get_ymin()-HopperLength)/2.0, HopperHeight);
                                B = Vec3D(0.0, (get_ymax()-get_ymin()-HopperExitLength)/2.0, HopperCornerHeight);
                                C = Vec3D(0.0, (get_ymax()-get_ymin()-HopperExitLength)/2.0, 0.0);
                
                
                
                                //now rotate the coordinates of A,B,C to be in (x,y,z) direction                
                                A = Vec3D(c*A.X-s*A.Z, A.Y, s*A.X+c*A.Z);
                                B = Vec3D(c*B.X-s*B.Z, B.Y, s*B.X+c*B.Z);
                                C = Vec3D(c*C.X-s*C.Z, C.Y, s*C.X+c*C.Z);
                                // the position of A determines shift and zmax
                                A.X +=shift;
                                B.X +=shift;
                                C.X +=shift;
                
                                //This lifts the hopper a distance above the chute
                                A.Z+=lift;
                                B.Z+=lift;
                                C.Z+=lift;
                                
                                
                                
                                //create a finite wall from B to A and from C to B
                                temp = B-A;
                                normal=Cross(Vec3D(-c,0,-s),temp)/sqrt(temp.GetLength2());
                                //normal  = Vec3D(0.0,temp.Z,-temp.Y) / sqrt(temp.GetLength2());
                                Walls[n+2].add_finite_wall(normal, Dot(normal,A));
                                temp = C-B;
                                //normal  = Vec3D(0.0,temp.Z,-temp.Y) / sqrt(temp.GetLength2());
                                normal=Cross(Vec3D(-c,0,-s),temp)/sqrt(temp.GetLength2());
                                Walls[n+2].add_finite_wall(normal, Dot(normal,B));
                                temp = A-C;
                                //normal = Vec3D(0.0,temp.Z,-temp.Y)/sqrt(temp.GetLength2());
                                normal=Cross(Vec3D(-c,0,-s),temp)/sqrt(temp.GetLength2());
                                Walls[n+2].add_finite_wall(normal,Dot(normal,C));
                                
                                
                                //Now for the right y-wall
                                A = Vec3D(0.0, (get_ymax()-get_ymin()+HopperLength)/2.0,HopperHeight);
                                B = Vec3D(0.0, (get_ymax()-get_ymin()+HopperExitLength)/2.0,HopperCornerHeight);
                                C = Vec3D(0.0, (get_ymax()-get_ymin()+HopperExitLength)/2.0,0.0);
                                
                                //now rotate the coordinates of A,B,C to be in (x,y,z) direction                
                                A = Vec3D(c*A.X-s*A.Z, A.Y, s*A.X+c*A.Z);
                                B = Vec3D(c*B.X-s*B.Z, B.Y, s*B.X+c*B.Z);
                                C = Vec3D(c*C.X-s*C.Z, C.Y, s*C.X+c*C.Z);
                                // the position of A determines shift and zmax
                                A.X +=shift;
                                B.X +=shift;
                                C.X +=shift;
                
                                //This lifts the hopper a distance above the chute
                                A.Z+=lift;
                                B.Z+=lift;
                                C.Z+=lift;
                                
                                //create a finite wall from B to A and from C to B
                                temp = A-B;
                                normal=Cross(Vec3D(-c,0,-s),temp)/sqrt(temp.GetLength2());
                                //normal  = Vec3D(0.0,-temp.Z,temp.Y) / sqrt(temp.GetLength2());
                                Walls[n+3].add_finite_wall(normal, Dot(normal,A));
                                temp = B-C;
                                //normal  = Vec3D(0.0,-temp.Z,temp.Y) / sqrt(temp.GetLength2());
                                normal=Cross(Vec3D(-c,0,-s),temp)/sqrt(temp.GetLength2());
                                Walls[n+3].add_finite_wall(normal, Dot(normal,B));
                                temp = C-A;
                                //normal = Vec3D(0.0,-temp.Z,temp.Y)/sqrt(temp.GetLength2());
                                normal=Cross(Vec3D(-c,0,-s),temp)/sqrt(temp.GetLength2());
                                Walls[n+3].add_finite_wall(normal,Dot(normal,C));
                                
                                
                                
                                
                                
                                
                        
                                
                        }
                
                
                
                
                
                
                
                //now shift the fixed particles at the bottom so that they begin where the chute begins
                for (vector<Particle*>::iterator it= get_ParticleHandler().begin(); it!=get_ParticleHandler().end(); ++it) { 
                        (*it)->get_Position().X+=shift;
                        #ifdef USE_SIMPLE_VERLET_INTEGRATION
                                (*it)->PrevPosition.X +=shift;
                        #endif
                }
        }
void ChuteWithHopper::constructor ( ) [inline]

This is the actually constructor, get called by all constructors above.

Reimplemented from Chute.

References align_base, centerHopper, fill_percent, hopper_dim, lift, lowerFillHeight, set_Hopper(), and shift.

Referenced by ChuteWithHopper().

                {
                lowerFillHeight=0.5;
                lift=0.0;
                set_Hopper(0.01, 0.01, 60.0, 0.08);
                shift = 0.0;    
                hopper_dim=1;
                align_base=true;

                fill_percent=50.0;
                centerHopper=false;

                }
virtual void ChuteWithHopper::create_inflow_particle ( ) [inline, virtual]

This creates an inflow particle in the top 50% of the hopper i.e.

between gamma=0.5 and gamma=1.0 Gamma is random number in the z direction and delta in the y direction In the 2D (hopper) case the particles are generated with equal probability in the y-direction, i.e. delta is from the edge of the domain In the 3D (hopper) case a third vector AD is generated and delta is again created for the sloping walls of the hopper

hopper_add_particle.jpg
Image shows the vectors in 2-dimension used to find a position inside the hopper
Bug:
for periodic walls this should be only minus one particle radius, this should be fixed at some point. Thomas' response: using one particle radius gives problems when the wall is not orthogonal to the y-direction; the distance has to be slightly higher than one; if you can figure out the exact value, then correct it please.

Reimplemented from Chute.

References centerHopper, Particle::compute_particle_mass(), fill_percent, mathsFunc::gamma(), Chute::get_ChuteAngle(), Particle::get_Position(), Particle::get_Radius(), RNG::get_RN(), MD::get_ymax(), MD::get_ymin(), hopper_dim, HopperAngle, HopperExitLength, HopperHeight, HopperLength, lift, Chute::MaxInflowParticleRadius, Chute::MinInflowParticleRadius, Chute::P0, MD::random, Particle::set_Position(), Particle::set_Radius(), MD::Species, and Vec3D::Z.

        {
                
                
                //use this formula to obtain bidispersed particles
                //P0.Radius = random(0.0,1.0)<0.1?MinInflowParticleRadius:MaxInflowParticleRadius;
                
                //the following formula yields polydispersed particle radii:

                P0.set_Radius(random.get_RN(MinInflowParticleRadius,MaxInflowParticleRadius));
                P0.compute_particle_mass(Species);
                
                //Define a orthogonal coordinate system this is usful in the hopper, see diagram in html documentation for details.
                static Mdouble s = sin(get_ChuteAngle());
                static Mdouble c = cos(get_ChuteAngle());
                static Mdouble Ht = tan(HopperAngle);
                static Mdouble Hc = cos(HopperAngle);
                static Vec3D AB = Vec3D(c,0.0,s);
                static Vec3D AC = Vec3D(-s,0.0,c);
                static Vec3D AD = Vec3D(0.0,1.0,0.0);
                
                //Point A is located in the centre of the hopper.
                static Vec3D A = Vec3D
                        (
                        centerHopper?40:0.0, 
                        (get_ymax()-get_ymin())/2.0, 
                        s*(-0.5*(HopperLength-HopperExitLength)) + c*HopperHeight
                        ) 
                        + AB*0.5*HopperLength 
                        + AC*(-0.5*HopperLength/Ht);

                Mdouble gamma = random.get_RN((100.0-fill_percent)/100.0,1.0);

                
                //              Mdouble gamma = random(lowerFillHeight,1.0);

                Mdouble delta;
                
                if (hopper_dim==1)
                        {
                        
                                
                                //For the one dimensional delta is a random distance between the two walls the -minus 2 particle radii is to stop 
                                //delta = random(ymin+P0.Radius,ymax-P0.Radius);

                                delta = random.get_RN(-0.5,0.5)*(get_ymax()-get_ymin()-2.0*P0.get_Radius());


                        }
                else
                        {
                                
                                delta= (random.get_RN(-1.0,1.0)*(0.5*gamma*HopperLength -P0.get_Radius()/Hc));

                        }
                P0.set_Position( A 
                        + AC * (gamma*0.5*HopperLength/Ht)

                        + AB * (random.get_RN(-1.0,1.0)*(0.5*gamma*HopperLength - P0.get_Radius()/Hc)) 
                        + AD*delta);
                        
                
                P0.get_Position().Z+=lift;
                
                //P0.Position.Y = random(ymin+P0.Radius, ymax-P0.Radius);
                //P0.Position.Y=delta;

                
                
        }

Allows chute length to be accessed.

Reimplemented from Chute.

References MD::get_xmax(), and shift.

{return get_xmax()-shift;}

References lift.

{return lift;}

Allows chute length to be accessed.

References Chute::ChuteAngle, MD::get_xmax(), HopperHeight, and shift.

                                                     {
          Mdouble height = HopperHeight+(get_xmax()-shift)*sin(ChuteAngle); 

          return sqrt(2.0*9.8*height);
        }
void ChuteWithHopper::lift_hopper ( Mdouble  distance) [inline]

This lifts the hopper above the plane of the chute.

References lift.

{lift=distance;}
virtual void ChuteWithHopper::print ( std::ostream &  os) [inline, virtual]

References HopperAngle, HopperExitHeight, HopperExitLength, HopperHeight, and HopperLength.

                                           {
                Chute::print(os); 
                os 
                  << ", HopperExitLength:" << HopperExitLength 
                        << ", HopperExitHeight:" << HopperExitHeight 
                        << ", HopperLength:" << HopperLength 
                        << ", HopperAngle:" << HopperAngle 
                        << ", HopperHeight:" << HopperHeight 
                        << endl;
        }
virtual void ChuteWithHopper::read ( std::istream &  is) [inline, virtual]

This function reads all chute data.

Reimplemented from Chute.

References HopperAngle, HopperExitHeight, HopperExitLength, HopperHeight, HopperLength, and shift.

int ChuteWithHopper::readNextArgument ( unsigned int &  i,
unsigned int  argc,
char *  argv[] 
) [inline, virtual]

Reimplemented from Chute.

References align_base, centerHopper, HopperAngle, HopperExitHeight, HopperExitLength, HopperHeight, HopperLength, lift, lowerFillHeight, and shift.

                                                                               {
                if (!strcmp(argv[i],"-hopperLength")) {
                        HopperLength=(atof(argv[i+1]));
                } else if (!strcmp(argv[i],"-hopperHeight")) {
                        HopperHeight=(atof(argv[i+1]));
                } else if (!strcmp(argv[i],"-hopperAngle")) {
                        HopperAngle=(atof(argv[i+1]));
                } else if (!strcmp(argv[i],"-hopperExitLength")) {
                        HopperExitLength=(atof(argv[i+1]));
                } else if (!strcmp(argv[i],"-hopperExitHeight")) {
                        HopperExitHeight=(atof(argv[i+1]));
                } else if (!strcmp(argv[i],"-lowerFillHeight")) {
                        lowerFillHeight=(atof(argv[i+1]));
                } else if (!strcmp(argv[i],"-centerHopper")) {
                        centerHopper=(atoi(argv[i+1]));
                } else if (!strcmp(argv[i],"-alignBase")) {
                        align_base=(atoi(argv[i+1]));
                } else if (!strcmp(argv[i],"-shift")) {
                        shift=(atof(argv[i+1]));
                } else if (!strcmp(argv[i],"-lift")) {
                        lift=(atof(argv[i+1]));
                } else return Chute::readNextArgument(i, argc, argv); //if argv[i] is not found, check the commands in Chute
                return true; //returns true if argv[i] is found
        }
void ChuteWithHopper::set_align_base ( bool  new_align) [inline]

References align_base.

{align_base=new_align;}
void ChuteWithHopper::set_centerHopper ( bool  new_) [inline]

References centerHopper.

{centerHopper=new_; }
void ChuteWithHopper::set_ChuteLength ( Mdouble  new_) [inline, virtual]

Allows chute length to be changed.

Reimplemented from Chute.

References MD::set_xmax(), MD::set_xmin(), and shift.

{if (new_>=0.0) {set_xmax(new_+shift); set_xmin(0.0);} else cerr << "WARNING : Chute length unchanged, value must be greater than or equal to zero" << endl;}
void ChuteWithHopper::set_Hopper ( Mdouble  ExitLength,
Mdouble  ExitHeight,
Mdouble  Angle,
Mdouble  Length 
) [inline]

References HopperAngle, HopperExitHeight, HopperExitLength, HopperLength, and constants::pi.

Referenced by constructor().

                                                                                              {
                if (ExitLength>=0.0) {HopperExitLength = ExitLength;} else cerr << "WARNING : Hopper exit length must be greater than or equal to zero" << endl;
                if (ExitHeight>=0.0) {HopperExitHeight = ExitHeight;} else cerr << "WARNING : Hopper exit height must be greater than or equal to zero" << endl;
                if (Angle>0.0&&Angle<90.0) {HopperAngle = Angle*constants::pi/180.0;} else cerr << "WARNING : Hopper angle must in (0,90)" << endl;
                if (Length>ExitLength) {HopperLength = Length;} else cerr << "WARNING : Hopper length must be greater than exit length" << endl;
        }
void ChuteWithHopper::set_hopper_dim ( Mdouble  new_hopper_dim) [inline]

References hopper_dim.

{hopper_dim=new_hopper_dim;}

References fill_percent.

{fill_percent=new_fill;}

References lowerFillHeight.

{lowerFillHeight=new_; }
void ChuteWithHopper::set_shift ( Mdouble  new_) [inline]

References MD::get_xmax(), MD::set_xmax(), and shift.

Referenced by add_hopper().

{if (new_>=0.0) {set_xmax(get_xmax()+new_-shift); shift = new_;} else cerr << "WARNING : Shift length unchanged, value must be greater than or equal to zero" << endl;}
virtual void ChuteWithHopper::setup_particles_initial_conditions ( ) [inline, virtual]

initialize particle position, velocity, radius

This initially set up the particles///.

Reimplemented from Chute.

References add_hopper().

        {
                Chute::setup_particles_initial_conditions();
                //cout << shift << " " << get_xmax() << " " << get_ParticleHandler().get_NumberOfParticles() << endl;
                add_hopper();   
        }       
virtual void ChuteWithHopper::write ( std::ostream &  os) [inline, virtual]

This function writes all chute data.

Reimplemented from Chute.

References HopperAngle, HopperExitHeight, HopperExitLength, HopperHeight, HopperLength, and shift.

                                            { 
                Chute::write(os); 
                os << HopperExitLength << " " << HopperExitHeight << " " << HopperLength 
                        << " " << HopperAngle << " " << HopperHeight << " " << shift << " " << endl;
        }

Member Data Documentation

This is the flag, which sets if the chute bottom is aligned with the hopper, by default it is.

Referenced by constructor(), readNextArgument(), and set_align_base().

If theis flag is set, the hopper will be constructed in the xy-center of the domain, and not next to the xmin-domain boundary; by default off.

Referenced by add_hopper(), constructor(), create_inflow_particle(), readNextArgument(), and set_centerHopper().

This is which percentage of the hopper is used for creating new partices;.

Referenced by constructor(), create_inflow_particle(), and set_HopperFillPercentage().

unsigned int ChuteWithHopper::hopper_dim [private]

This is the dimension of the hopper, my default it is one dimensional and hence does not have side wall.

Referenced by add_hopper(), constructor(), create_inflow_particle(), and set_hopper_dim().

Angle between the two pieces of the hopper walls.

Referenced by add_hopper(), create_inflow_particle(), print(), read(), readNextArgument(), set_Hopper(), and write().

Dimension of the hopper exit in vertical direction.

Referenced by add_hopper(), print(), read(), readNextArgument(), set_Hopper(), and write().

Dimension of the hopper exit in vertical direction.

Referenced by add_hopper(), create_inflow_particle(), print(), read(), readNextArgument(), set_Hopper(), and write().

Dimension of the hopper in horizontal direction.

Referenced by add_hopper(), create_inflow_particle(), get_MaximumVelocityInducedByGravity(), print(), read(), readNextArgument(), and write().

Dimension of the hopper in vertical direction.

Referenced by add_hopper(), create_inflow_particle(), print(), read(), readNextArgument(), set_Hopper(), and write().

This is the vertical distance the chute is lifted above the plane.

Referenced by add_hopper(), constructor(), create_inflow_particle(), get_lift_hopper(), lift_hopper(), and readNextArgument().

Relative height (in [0,1)) above which teh hopper is replenished with new particles.

Referenced by constructor(), readNextArgument(), and set_lowerFillHeight().

The x position where the Chute starts (defined as the beginning of the hopper)

Referenced by add_hopper(), constructor(), get_ChuteLength(), get_MaximumVelocityInducedByGravity(), read(), readNextArgument(), set_ChuteLength(), set_shift(), and write().


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