MercuryDPM  0.10
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
HGRID_3D Class Reference

This adds on the hierarchical grid code for 3D problems. More...

#include <HGRID_3D.h>

+ Inheritance diagram for HGRID_3D:

Public Member Functions

 HGRID_3D ()
 This is the default constructor. All it does is set senible defaults. More...
 
 HGRID_3D (MD &other)
 Copy-constructor for creates an HGRID problem from an existing MD problem. More...
 
 HGRID_3D (HGRID_base &other)
 
void constructor ()
 This is the actually constructor it is called do both constructors above. More...
 
- Public Member Functions inherited from HGRID_base
 HGRID_base ()
 This is the default constructor. All it does is set senible defaults. More...
 
 ~HGRID_base ()
 This is the default destructor. More...
 
 HGRID_base (MD &other)
 Copy-constructor for creates an HGRID problem from an existing MD problem. More...
 
void constructor ()
 This is the actually constructor it is called do both constructors above. More...
 
void HGRID_actions_before_time_loop ()
 This sets up the broad phase information, has to be done at this stage becuase it requires the partcle size. More...
 
void HGRID_actions_before_time_step ()
 This resets all the bucket information. More...
 
void set_HGRID_num_buckets (unsigned int new_num_buckets)
 This sets the number of buckets for the HGRID. More...
 
void set_HGRID_num_buckets_to_power ()
 set number of buckets to the smallest power of two bigger than the number of particles More...
 
void set_HGRID_num_buckets_to_power (unsigned int N)
 set number of buckets to the smallest power of two bigger than N More...
 
void read (std::istream &is)
 This function reads all HGRID data. More...
 
void write (std::ostream &os)
 This function writes all HGRID data. More...
 
void print (std::ostream &os, bool print_all)
 This function outputs all HGRID data. More...
 
Mdouble getHGridCurrentMaxRelativeDisplacement ()
 
Mdouble getHGridTotalCurrentMaxRelativeDisplacement ()
 
void setHGridUpdateEachTimeStep (bool updateEachTimeStep)
 
bool getHGridUpdateEachTimeStep ()
 
void setHGridMaxLevels (int HGridMaxLevels)
 
int getHGridMaxLevels ()
 
HGridMethod getHGridMethod ()
 
void setHGridMethod (HGridMethod hGridMethod)
 
HGridDistribution getHGridDistribution ()
 
void setHGridDistribution (HGridDistribution hGridDistribution)
 
Mdouble getHGridCellOverSizeRatio ()
 
void setHGridCellOverSizeRatio (Mdouble cellOverSizeRatio)
 
- Public Member Functions inherited from MD
void constructor ()
 A public constructor, which sets defaults so the problem can be solved off the shelf. More...
 
 MD ()
 
 MD (STD_save &other)
 
virtual ~MD ()
 
void info ()
 Set up a virtual info this will be provided from the inhertiance. More...
 
void solve ()
 The work horse of the code. More...
 
void solve (unsigned int argc, char *argv[])
 Read arguments before solving. More...
 
void solveWithMDCLR ()
 Tries to solve the problem using MDCLR. More...
 
Mdouble get_t ()
 Access function for the time. More...
 
void set_t (Mdouble new_t)
 Access function for the time. More...
 
int get_NSpecies ()
 Allows the number of Species to be accessed. More...
 
std::vector< CSpecies > & get_Species (void)
 Allows the species to be copied. More...
 
CSpeciesget_Species (int i)
 Allows the species to be accessed. More...
 
CSpeciesget_MixedSpecies (int i, int j)
 Allows the mixed species to be accessed. More...
 
void set_MixedSpecies (int i, int j, CSpecies &S)
 Allows the mixed species to be set. More...
 
void set_tmax (Mdouble new_tmax)
 Allows the upper time limit to be changed. More...
 
Mdouble get_tmax ()
 Allows the upper time limit to be accessed. More...
 
ParticleHandlergetParticleHandler ()
 
WallHandlergetWallHandler ()
 
BoundaryHandlergetBoundaryHandler ()
 
void set_savecount (int new_)
 Allows the number of time steps between saves to be changed, see also set_number_of_saves. More...
 
void set_save_count_all (int new_)
 
void set_save_count_data (int new_)
 
void set_save_count_ene (int new_)
 
void set_save_count_stat (int new_)
 
void set_save_count_fstat (int new_)
 
int get_savecount ()
 Allows the number of time steps between saves to be accessed. More...
 
int get_save_count ()
 
int get_save_count_data ()
 
int get_save_count_ene ()
 
int get_save_count_stat ()
 
int get_save_count_fstat ()
 
void set_do_stat_always (bool new_)
 Sets how often the data is saved using the number of saves wanted, tmax, and dt. See also set_savecount. More...
 
void set_number_of_saves (Mdouble N)
 
void set_number_of_saves_all (Mdouble N)
 
void set_number_of_saves_data (Mdouble N)
 
void set_number_of_saves_ene (Mdouble N)
 
void set_number_of_saves_stat (Mdouble N)
 
void set_number_of_saves_fstat (Mdouble N)
 
void set_plastic_k1_k2max_kc_depth (Mdouble k1_, Mdouble k2max_, Mdouble kc_, Mdouble depth_, unsigned int indSpecies=0)
 Allows the plastic constants to be changed. More...
 
void set_k1 (Mdouble new_, unsigned int indSpecies=0)
 
void set_k2max (Mdouble new_, unsigned int indSpecies=0)
 
void set_kc (Mdouble new_, unsigned int indSpecies=0)
 
void set_depth (Mdouble new_, unsigned int indSpecies=0)
 
Mdouble get_k1 (unsigned int indSpecies=0)
 Allows the plastic constants to be accessed. More...
 
Mdouble get_k2max (unsigned int indSpecies=0)
 
Mdouble get_kc (unsigned int indSpecies=0)
 
Mdouble get_depth (unsigned int indSpecies=0)
 
Mdouble get_plastic_dt (Mdouble mass, unsigned int indSpecies=0)
 
void set_k (Mdouble new_, unsigned int indSpecies=0)
 Allows the spring constant to be changed. More...
 
Mdouble get_k (int indSpecies=0)
 Allows the spring constant to be accessed. More...
 
void set_kt (Mdouble new_, unsigned int indSpecies=0)
 Allows the spring constant to be changed. More...
 
Mdouble get_kt (int indSpecies=0)
 Allows the spring constant to be accessed. More...
 
void set_krolling (Mdouble new_, unsigned int indSpecies=0)
 Allows the spring constant to be changed. More...
 
Mdouble get_krolling (int indSpecies=0)
 Allows the spring constant to be accessed. More...
 
void set_ktorsion (Mdouble new_, unsigned int indSpecies=0)
 Allows the spring constant to be changed. More...
 
Mdouble get_ktorsion (int indSpecies=0)
 Allows the spring constant to be accessed. More...
 
void set_rho (Mdouble new_, unsigned int indSpecies=0)
 Allows the density to be changed. More...
 
Mdouble get_rho (int indSpecies=0)
 Allows the density to be accessed. More...
 
void set_dispt (Mdouble new_, unsigned int indSpecies=0)
 Allows the tangential viscosity to be changed. More...
 
Mdouble get_dispt (unsigned int indSpecies=0)
 Allows the tangential viscosity to be accessed. More...
 
void set_disprolling (Mdouble new_, unsigned int indSpecies=0)
 Allows the tangential viscosity to be changed. More...
 
Mdouble get_disprolling (unsigned int indSpecies=0)
 Allows the tangential viscosity to be accessed. More...
 
void set_disptorsion (Mdouble new_, unsigned int indSpecies=0)
 Allows the tangential viscosity to be changed. More...
 
Mdouble get_disptorsion (unsigned int indSpecies=0)
 Allows the tangential viscosity to be accessed. More...
 
void set_disp (Mdouble new_, unsigned int indSpecies=0)
 Allows the normal dissipation to be changed. More...
 
Mdouble get_disp (unsigned int indSpecies=0)
 Allows the normal dissipation to be accessed. More...
 
void set_dissipation (Mdouble new_, unsigned int indSpecies=0)
 Allows the normal dissipation to be changed. More...
 
Mdouble get_dissipation (unsigned int indSpecies=0)
 Allows the normal dissipation to be accessed. More...
 
void set_mu (Mdouble new_, unsigned int indSpecies=0)
 Allows the Coulomb friction coefficient to be changed. More...
 
Mdouble get_mu (unsigned int indSpecies=0)
 Allows the Coulomb friction coefficient to be accessed. More...
 
void set_murolling (Mdouble new_, unsigned int indSpecies=0)
 Allows the Coulomb friction coefficient to be changed. More...
 
Mdouble get_murolling (unsigned int indSpecies=0)
 Allows the Coulomb friction coefficient to be accessed. More...
 
void set_mutorsion (Mdouble new_, unsigned int indSpecies=0)
 Allows the Coulomb friction coefficient to be changed. More...
 
Mdouble get_mutorsion (unsigned int indSpecies=0)
 Allows the Coulomb friction coefficient to be accessed. More...
 
void set_rotation (bool new_)
 
bool get_rotation ()
 
void set_dim_particle (int new_, unsigned int indSpecies=0)
 Allows the dimension of the particle (f.e. for mass) to be changed. More...
 
int get_dim_particle (unsigned int indSpecies=0)
 Allows the dimension of the particle (f.e. for mass) to be accessed. More...
 
bool get_save_data_data ()
 Returns the data counter. More...
 
bool get_save_data_ene ()
 
bool get_save_data_fstat ()
 
bool get_save_data_stat ()
 
bool get_do_stat_always ()
 
void set_k_and_restitution_coefficient (Mdouble k_, Mdouble eps, Mdouble mass, unsigned int indSpecies=0)
 Sets k, disp such that it matches a given tc and eps for a collision of two copies of P. More...
 
void set_collision_time_and_restitution_coefficient (Mdouble tc, Mdouble eps, Mdouble mass, unsigned int indSpecies=0)
 Sets k, disp such that it matches a given tc and eps for a collision of two copies of P. More...
 
void set_collision_time_and_restitution_coefficient (Mdouble tc, Mdouble eps, Mdouble mass1, Mdouble mass2, unsigned int indSpecies=0)
 Set k, disp such that is matches a given tc and eps for a collision of two different masses. More...
 
void set_collision_time_and_normal_and_tangential_restitution_coefficient (Mdouble tc, Mdouble eps, Mdouble beta, Mdouble mass1, Mdouble mass2, unsigned int indSpecies=0)
 See CSpecies::set_collision_time_and_normal_and_tangential_restitution_coefficient. More...
 
void set_collision_time_and_normal_and_tangential_restitution_coefficient_nodispt (Mdouble tc, Mdouble eps, Mdouble beta, Mdouble mass1, Mdouble mass2, unsigned int indSpecies=0)
 See CSpecies::set_collision_time_and_normal_and_tangential_restitution_coefficient. More...
 
Mdouble get_collision_time (Mdouble mass, unsigned int indSpecies=0)
 Calculates collision time for two copies of a particle of given disp, k, mass. More...
 
Mdouble get_restitution_coefficient (Mdouble mass, unsigned int indSpecies=0)
 Calculates restitution coefficient for two copies of given disp, k, mass. More...
 
Mdouble get_xmin ()
 Get xmin. More...
 
Mdouble get_xmax ()
 Get xmax. More...
 
Mdouble get_ymin ()
 Gets ymin. More...
 
Mdouble get_ymax ()
 Gets ymax. More...
 
Mdouble get_zmin ()
 Gets zmin. More...
 
Mdouble get_zmax ()
 Gets zmax. More...
 
void set_xmin (Mdouble new_xmin)
 Sets xmin and walls, assuming the standard definition of walls as in the default constructor. More...
 
void set_ymin (Mdouble new_ymin)
 
void set_zmin (Mdouble new_zmin)
 Sets ymin and walls, assuming the standard definition of walls as in the default constructor. More...
 
void set_xmax (Mdouble new_xmax)
 Sets xmax and walls, assuming the standard definition of walls as in the default constructor. More...
 
void set_ymax (Mdouble new_ymax)
 Sets ymax and walls, assuming the standard definition of walls as in the default constructor. More...
 
void set_zmax (Mdouble new_zmax)
 Sets ymax and walls, assuming the standard definition of walls as in the default constructor. More...
 
void set_dt (Mdouble new_dt)
 Allows the time step dt to be changed. More...
 
Mdouble get_dt ()
 Allows the time step dt to be accessed. More...
 
void set_name (const char *name)
 Sets the name of the problem, used for the same data files. More...
 
void set_xballs_colour_mode (int new_cmode)
 Set the xball output mode. More...
 
void set_xballs_cmode (int new_cmode)
 
int get_xballs_cmode ()
 
void set_xballs_vector_scale (double new_vscale)
 Set the scale of vectors in xballs. More...
 
double get_xballs_vscale ()
 
void set_xballs_additional_arguments (std::string new_)
 Set the additional arguments for xballs. More...
 
std::string get_xballs_additional_arguments ()
 
void set_xballs_scale (Mdouble new_scale)
 Set the scale of the xballs problem. The default is fit to screen. More...
 
double get_xballs_scale ()
 
void set_gravity (Vec3D new_gravity)
 Allows the gravitational acceleration to be changed. More...
 
Vec3D get_gravity ()
 Allows the gravitational acceleration to be accessed. More...
 
void set_dim (int new_dim)
 Allows the dimension of the simulation to be changed. More...
 
int get_dim ()
 Allows the dimension of the simulation to be accessed. More...
 
int get_restart_version ()
 Gets restart_version. More...
 
void set_restart_version (int new_)
 Sets restart_version. More...
 
bool get_restarted ()
 Gets restarted. More...
 
Mdouble get_max_radius ()
 Sets restarted. More...
 
void set_restarted (bool new_)
 
bool get_append ()
 Gets restarted. More...
 
void set_append (bool new_)
 Sets restarted. More...
 
Mdouble get_ene_ela ()
 Gets ene_ela. More...
 
void set_ene_ela (Mdouble new_)
 Sets ene_ela. More...
 
void add_ene_ela (Mdouble new_)
 Sets ene_ela. More...
 
void Remove_Particle (int IP)
 This function removes partice IP from the vector of particles by moving the last particle in the vector to the position if IP Also it checks if the moved Particle has any tangentialsspring-information, which needs to be moved to a different particle, because tangential spring information always needs to be stored in the real particle with highest particle index. More...
 
Mdouble get_Mass_from_Radius (Mdouble radius, int indSpecies=0)
 
Mdouble get_maximum_velocity (BaseParticle &P)
 Calculates the maximum velocity allowed for a collision of two copies of P (for higher velocities particles could pass through each other) More...
 
virtual BaseParticlegetSmallestParticle ()
 
virtual BaseParticlegetLargestParticle ()
 
virtual void removeParticle (int iP)
 
Mdouble get_maximum_velocity ()
 
void set_dt_by_mass (Mdouble mass)
 Sets dt to 1/50-th of the collision time for two particles of mass P. More...
 
void set_dt (BaseParticle &P)
 Sets dt to 1/50-th of the collision time for two copies of P. More...
 
void set_dt ()
 Sets dt to 1/50-th of the smallest possible collision time. More...
 
virtual void setup_particles_initial_conditions ()
 This function allows the initial conditions of the particles to be set, by default locations is random. More...
 
virtual void create_xballs_script ()
 This creates a scipt which can be used to load the xballs problem to display the data just generated. More...
 
virtual double getInfo (BaseParticle &P)
 Allows the user to set what is written into the info column in the data file. By default is store the Species ID number. More...
 
virtual void save_restart_data ()
 Stores all MD data. More...
 
int load_restart_data ()
 Loads all MD data. More...
 
int load_restart_data (std::string filename)
 
void statistics_from_restart_data (const char *name)
 Loads all MD data and plots statistics for all timesteps in the .data file. More...
 
virtual void write_v1 (std::ostream &os)
 Writes all MD data. More...
 
virtual void read_v1 (std::istream &is)
 Reads all MD data. More...
 
virtual void read_v2 (std::istream &is)
 
bool load_from_data_file (const char *filename, unsigned int format=0)
 This allows particle data to be reloaded from data files. More...
 
bool load_par_ini_file (const char *filename)
 allows the user to read par.ini files (useful to read MDCLR files) More...
 
bool read_next_from_data_file (unsigned int format=0)
 by default format do not pass an argument; only specify format if you have to read a special format (f.e. dim=2, but format=14 (3d format)) More...
 
int read_dim_from_data_file ()
 
bool find_next_data_file (Mdouble tmin, bool verbose=true)
 
void add_Species (CSpecies &S)
 
void add_Species (void)
 
void set_format (int new_)
 
int get_format ()
 
int readArguments (unsigned int argc, char *argv[])
 Can interpret main function input arguments that are passed by the driver codes. More...
 
- Public Member Functions inherited from STD_save
 STD_save ()
 Default constructor: sets the counter to 0 (i.e. no number will be included). More...
 
 STD_save (STD_save &other)
 Copy constructor. More...
 
void constructor ()
 
void inc_counter_in_file ()
 Increament the counter value stored in the file_counter by 1 and store the new value. More...
 
int read_run_num_from_file ()
 Read rom the counter file the counter. More...
 
void set_counter_from_file ()
 Sets the counter based on the current number stored in the counter file. More...
 
void save_info_to_disk ()
 Saves the information generated by info to disk in a file. More...
 
void set_counter (int new_counter)
 This set the counter, overriding the defaults. More...
 
int get_counter ()
 This returns the current value of the counter. More...
 
bool FileExists (std::string strFilename)
 Function to check if a file exists, is used to check if a run has already need done. More...
 
void auto_number ()
 
std::vector< intget_numbers (int size_x, int size_y)
 This turns a counter into two indexs for doing parmater studies. The indexs run from 1:size_x and 1:size_y where as the study number starts at 0. More...
 
int launch_new (const char *name, bool quick=false)
 This launch a code from within this code. Please pass the name of the code to run. More...
 
void set_name (const char *name)
 Sets the name of the problem, used for the same data files. More...
 
std::string get_name ()
 Allows the problem_name to be accessed. More...
 
std::fstream & get_data_file ()
 Allows the problem_name to be accessed. More...
 
std::fstream & get_stat_file ()
 Allows the problem_name to be accessed. More...
 
std::fstream & get_fstat_file ()
 Allows the problem_name to be accessed. More...
 
std::fstream & get_ene_file ()
 Allows the problem_name to be accessed. More...
 
void set_fstat_filename (std::string filename)
 
void set_data_filename (std::string filename)
 
void set_stat_filename (std::string filename)
 
void set_ene_filename (std::string filename)
 
void set_fstat_filename ()
 
void set_data_filename ()
 
void set_stat_filename ()
 
void set_ene_filename ()
 
std::string get_fstat_filename ()
 
std::string get_data_filename ()
 
std::string get_stat_filename ()
 
std::string get_ene_filename ()
 
void set_step_size (unsigned int new_)
 
unsigned int get_step_size ()
 
void set_options_fstat (unsigned int new_)
 set and get for file options More...
 
unsigned int get_options_fstat (void)
 
void set_options_data (unsigned int new_)
 
unsigned int get_options_data (void)
 
void set_options_stat (unsigned int new_)
 
unsigned int get_options_stat (void)
 
void set_options_restart (unsigned int new_)
 
unsigned int get_options_restart (void)
 
void set_options_ene (unsigned int new_)
 
unsigned int get_options_ene (void)
 
bool open_file (std::fstream &file, std::string filename, unsigned int options, std::fstream::openmode mode)
 
bool open_fstat_file (std::fstream::openmode mode=std::fstream::out)
 
bool open_data_file (std::fstream::openmode mode=std::fstream::out)
 
bool open_stat_file (std::fstream::openmode mode=std::fstream::out)
 
bool open_ene_file (std::fstream::openmode mode=std::fstream::out)
 
bool open_counted_file (std::fstream &file, std::string filenameNoCount, std::fstream::openmode mode)
 opens file needed if data is written in multiple files More...
 
bool increase_counter_fstat (std::fstream::openmode mode)
 
bool increase_counter_data (std::fstream::openmode mode)
 
bool increase_counter_stat (std::fstream::openmode mode)
 
bool increase_counter_ene (std::fstream::openmode mode)
 
void set_file_counter (int new_)
 
int get_file_counter ()
 

Protected Member Functions

void HGRID_UpdateParticleInHgrid (BaseParticle *obj)
 This adds a partcile to the Grid, called in the grid setup routies. More...
 
void HGRID_RemoveParticleFromHgrid (BaseParticle *obj)
 
virtual void CheckCell (int x, int y, int z, int l, BaseParticle *obj, HGrid *grid)
 Check collisions for a general cell. More...
 
virtual void CheckCell_current (int x, int y, int z, int l, HGrid *grid)
 Checks for a collision in the particles own cell. More...
 
void CheckObjAgainstGrid (HGrid *grid, BaseParticle *obj)
 Check if an Particle has a collision in the grid; avoids multiple checks. More...
 
void CheckObjAgainstWholeGrid (HGrid *grid, BaseParticle *obj)
 Check if an Particle has a collision in the grid. More...
 
bool TestCell (int x, int y, int z, int l, BaseParticle *obj, HGrid *grid)
 Tests obj against all particles in cell similar to CheckCell, but links to TestObject instead of compute_internal_forces. More...
 
bool TestObjAgainstGrid (HGrid *grid, BaseParticle *obj)
 Tests obj against all neighbouring particles similar to CheckObjAgainstGrid, but links to TestCell instead of CheckCell. More...
 
- Protected Member Functions inherited from HGRID_base
void InitBroadPhase ()
 This sets up the parameters required for the contact model. More...
 
void HGRID_InsertParticleToHgrid (BaseParticle *obj)
 Inserts a single Particle to current grid. More...
 
void broad_phase (BaseParticle *i)
 This makes the board_phase of contact point at the HGRID code. More...
 
virtual bool TestObject (BaseParticle *pI, BaseParticle *pJ)
 criterium for inserting a particle (returns false, if particles overlap;) More...
 
void HGRID_update_move (BaseParticle *iP, Mdouble move)
 
void HGRID_actions_before_integration ()
 
void HGRID_actions_after_integration ()
 
int readNextArgument (unsigned int &i, unsigned int argc, char *argv[])
 
- Protected Member Functions inherited from MD
virtual void compute_all_forces ()
 This does the force computation. More...
 
virtual void compute_internal_forces (BaseParticle *i)
 Computes the forces between particles (internal in the sence that the sum over all these forces is zero i.e. fully modelled forces) More...
 
CTangentialSpringgetTangentialSpring (BaseParticle *PI, BaseParticle *PJ, BaseParticle *PJreal)
 
CTangentialSpringgetTangentialSpringWall (BaseParticle *pI, int w)
 
virtual void compute_internal_forces (BaseParticle *P1, BaseParticle *P2)
 Computes the forces between particles (internal in the sence that the sum over all these forces is zero i.e. fully modelled forces) More...
 
void compute_plastic_internal_forces (BaseParticle *P1, BaseParticle *P2)
 Computes plastic forces between particles. More...
 
virtual void compute_external_forces (BaseParticle *PI)
 This is were the computation of external forces takes place (e.g. gravity) More...
 
virtual void compute_walls (BaseParticle *PI)
 This is were the walls are. More...
 
Mdouble computeShortRangeForceWithWall (BaseParticle *pI, int wall, CSpecies *pSpecies, Mdouble dist)
 
Mdouble computeShortRangeForceWithParticle (BaseParticle *PI, BaseParticle *PJ, BaseParticle *PJreal, CSpecies *pSpecies, Mdouble dist)
 
virtual void actions_before_time_loop ()
 This is actions before the start of the main time loop. More...
 
virtual void HGRID_InsertParticleToHgrid (BaseParticle *obj UNUSED)
 This is action before the time step is started. More...
 
virtual void HGRID_UpdateParticleInHgrid (BaseParticle *obj UNUSED)
 
virtual void HGRID_RemoveParticleFromHgrid (BaseParticle *obj UNUSED)
 
virtual bool get_HGRID_UpdateEachTimeStep ()
 
virtual void actions_before_time_step ()
 This is action before the time step is started. More...
 
virtual void actions_after_solve ()
 This is actions at the end of the code, but before the files are closed. More...
 
virtual void actions_after_time_step ()
 This is action after the time step is finished. More...
 
virtual void output_xballs_data ()
 Output xball data for Particle i. More...
 
virtual void output_xballs_data_particle (int i)
 This function outputs the location and velocity of the particle in a format the xballs progream can read. More...
 
virtual void start_ene ()
 Functions for ene file. More...
 
virtual void fstat_header ()
 
virtual void output_ene ()
 This function outputs statistical data - The default is to compute the rotational kinetic engergy, linear kinetic energy, and the centre of mass. More...
 
virtual void initialize_statistics ()
 Functions for statistics. More...
 
virtual void output_statistics ()
 
virtual void gather_statistics_collision (int index1 UNUSED, int index2 UNUSED, Vec3D Contact UNUSED, Mdouble delta UNUSED, Mdouble ctheta UNUSED, Mdouble fdotn UNUSED, Mdouble fdott UNUSED, Vec3D P1_P2_normal_ UNUSED, Vec3D P1_P2_tangential UNUSED)
 
virtual void process_statistics (bool usethese UNUSED)
 
virtual void finish_statistics ()
 
virtual void set_initial_pressures_for_pressure_controlled_walls ()
 
virtual void do_integration_before_force_computation (BaseParticle *pI)
 This is were the integration is done. More...
 
virtual void checkInteractionWithBoundaries ()
 
virtual void do_integration_after_force_computation (BaseParticle *pI)
 This is were the integration is done. More...
 
void set_FixedParticles (unsigned int n)
 
void initialize_tangential_springs ()
 
void compute_particle_masses ()
 Computes the mass of each particle. More...
 
virtual void cout_time ()
 std::couts time More...
 
virtual bool continue_solve ()
 
void reset_DeltaMax ()
 sets the history parameter DeltaMax of all particles to zero More...
 
void reset_TangentialSprings ()
 sets the history parameter TangentialSprings of all particles to zero More...
 

Additional Inherited Members

- Public Attributes inherited from HGRID_base
HGridgrid
 
- Public Attributes inherited from MD
RNG random
 
- Protected Attributes inherited from MD
std::vector< CSpeciesSpecies
 These are the particle parameters like dissipation etc. More...
 
- Protected Attributes inherited from STD_save
std::stringstream problem_name
 Stores the problem_name. More...
 
std::stringstream data_filename
 These store the save file names, by default they are derived from problem_name. More...
 
std::stringstream stat_filename
 
std::stringstream fstat_filename
 
std::stringstream ene_filename
 
std::fstream data_file
 Stream used for data files. More...
 
std::fstream stat_file
 
std::fstream fstat_file
 
std::fstream ene_file
 
unsigned int options_fstat
 Indicators if files are created or not 0: file will not be created 1: file will be written in one file 2: file will be written in multiple files. More...
 
unsigned int options_data
 
unsigned int options_stat
 
unsigned int options_ene
 
unsigned int options_restart
 
unsigned int file_counter
 Counter needed if file will be written in multiple files. More...
 
unsigned int step_size
 

Detailed Description

This adds on the hierarchical grid code for 3D problems.

Definition at line 35 of file HGRID_3D.h.

Constructor & Destructor Documentation

HGRID_3D::HGRID_3D ( )
inline

This is the default constructor. All it does is set senible defaults.

Definition at line 40 of file HGRID_3D.h.

References constructor().

41  {
42  constructor();
43  #ifdef CONSTUCTOR_OUTPUT
44  std::cerr << "HGRID_3D() finished"<<std::endl;
45  #endif
46  }
void constructor()
This is the actually constructor it is called do both constructors above.
Definition: HGRID_3D.h:69
HGRID_3D::HGRID_3D ( MD other)
inline

Copy-constructor for creates an HGRID problem from an existing MD problem.

Bug:
Why is the constructor call commented out. This does not make sence to be Anthony. Dinant: I think this has something to do with somebody copying a 2D MD problem into a 3D Hgrid, then the 3D_Hgrid constructor changes properties of the MD class. Solution in my opinion is to change the HGRID_3D::constructor()

Definition at line 53 of file HGRID_3D.h.

53  : MD(other), HGRID_base(other)
54  {
55  /*constructor();*/
56  #ifdef CONSTUCTOR_OUTPUT
57  std::cerr << "HGRID_3D(MD& other) finished"<<std::endl;
58  #endif
59  }
MD()
Definition: MD.h:75
HGRID_base()
This is the default constructor. All it does is set senible defaults.
Definition: HGRID_base.cc:28
HGRID_3D::HGRID_3D ( HGRID_base other)
inline

Definition at line 60 of file HGRID_3D.h.

60  : MD(other), HGRID_base(other)
61  {
62  /*constructor();*/
63  #ifdef CONSTUCTOR_OUTPUT
64  std::cerr << "HGRID_3D(HGRID_base& other) finished"<<std::endl;
65  #endif
66  }
MD()
Definition: MD.h:75
HGRID_base()
This is the default constructor. All it does is set senible defaults.
Definition: HGRID_base.cc:28

Member Function Documentation

void HGRID_3D::CheckCell ( int  x,
int  y,
int  z,
int  l,
BaseParticle obj,
HGrid grid 
)
protectedvirtual

Check collisions for a general cell.

Check for collisions with a general cell.

Todo:
{What is the use of this check?}
Bug:
{TW: This check is not necessary, I believe. This is the most-expensive function in most codes (the two checks in this function slows down granular jet by 15%) and the selftests are not affected. DK: I do think this is neccesary, for example: If two cells hash to the same bucket and a particle in one of these cells check for collisions with the other cell. Then due to the hashingcollision it also gets all particles in it's own cell and thus generating false collisions.}

Definition at line 58 of file HGRID_3D.cc.

References MD::compute_internal_forces(), HGrid::ComputeHashBucketIndex(), BaseParticle::get_HGRID_Level(), BaseParticle::get_HGRID_NextObject(), BaseParticle::get_HGRID_x(), BaseParticle::get_HGRID_y(), BaseParticle::get_HGRID_z(), and HGrid::objectBucket.

Referenced by CheckObjAgainstGrid().

59 {
60  int bucket = grid->ComputeHashBucketIndex(x,y,z,l);
61 
63  if ((obj->get_HGRID_x() == x) && (obj->get_HGRID_y() == y) && (obj->get_HGRID_z() == z) && (obj->get_HGRID_Level() == l)) { return; }
64 
65  // Loop through all objects in the bucket to find nearby objects
66  BaseParticle *p = grid->objectBucket[bucket];
67 
68  while (p!=NULL)
69  {
71  //Check if the particle realy is in the target cell (i.e. no hashing error has occured)
72  if ((p->get_HGRID_x() == x) && (p->get_HGRID_y() == y) && (p->get_HGRID_z() == z) && (p->get_HGRID_Level() == l))
73  {
75  }
76  p = p->get_HGRID_NextObject();
77  }
78 }
int get_HGRID_z() const
BaseParticle * get_HGRID_NextObject() const
int get_HGRID_x() const
int ComputeHashBucketIndex(int x, int y, int z, int l)
Computes hash bucket index in range [0, NUM_BUCKETS-1].
Definition: HGRID.cc:93
int get_HGRID_y() const
int get_HGRID_Level() const
virtual void compute_internal_forces(BaseParticle *i)
Computes the forces between particles (internal in the sence that the sum over all these forces is ze...
Definition: MD.cc:1980
BaseParticle ** objectBucket
objectBucket[b] stores pointer to first element in hash bucket b; initially all NULL ...
Definition: HGRID.h:58
void HGRID_3D::CheckCell_current ( int  x,
int  y,
int  z,
int  l,
HGrid grid 
)
protectedvirtual

Checks for a collision in the particles own cell.

Checks for collision in the particles own cell.

Bug:
{TW: This check is not necessary, I believe. This is the most-expensive function in most codes (the two checks in this function slows down granular jet by 15%) and the selftests are not affected. DK: I do think this is neccesary, for example: If two cells hash to the same bucket and a particle in one of these cells check for collisions with the other cell. Then due to the hashingcollision it also gets all particles in it's own cell and thus generating false collisions.}

Definition at line 31 of file HGRID_3D.cc.

References HGrid::bucketIsChecked, MD::compute_internal_forces(), HGrid::ComputeHashBucketIndex(), BaseParticle::get_HGRID_Level(), BaseParticle::get_HGRID_NextObject(), BaseParticle::get_HGRID_x(), BaseParticle::get_HGRID_y(), BaseParticle::get_HGRID_z(), and HGrid::objectBucket.

Referenced by CheckObjAgainstGrid().

32 {
33  int bucket = grid->ComputeHashBucketIndex(x,y,z,l);
34  if (grid->bucketIsChecked[bucket]) { return; }
35  BaseParticle* p1 = grid->objectBucket[bucket];
36 
37  while (p1!=NULL)
38  {
40  while (p2!=NULL)
41  {
43  //Check if the particle realy is in the target cell (i.e. no hashing error has occured)
44  if ((p1->get_HGRID_x() == p2->get_HGRID_x()) && (p1->get_HGRID_y() == p2->get_HGRID_y()) && (p1->get_HGRID_z() == p2->get_HGRID_z()) && (p1->get_HGRID_Level() == p2->get_HGRID_Level()))
45  {
47  }
48  p2 = p2->get_HGRID_NextObject();
49  }
50  p1 = p1->get_HGRID_NextObject();
51  }
52  grid->bucketIsChecked[bucket] = true;
53 }
int get_HGRID_z() const
BaseParticle * get_HGRID_NextObject() const
int get_HGRID_x() const
int ComputeHashBucketIndex(int x, int y, int z, int l)
Computes hash bucket index in range [0, NUM_BUCKETS-1].
Definition: HGRID.cc:93
bool * bucketIsChecked
bucketIsChecked[b] stores if hash bucket b is checked already; initially all zero ...
Definition: HGRID.h:61
int get_HGRID_y() const
int get_HGRID_Level() const
virtual void compute_internal_forces(BaseParticle *i)
Computes the forces between particles (internal in the sence that the sum over all these forces is ze...
Definition: MD.cc:1980
BaseParticle ** objectBucket
objectBucket[b] stores pointer to first element in hash bucket b; initially all NULL ...
Definition: HGRID.h:58
void HGRID_3D::CheckObjAgainstGrid ( HGrid grid,
BaseParticle obj 
)
protectedvirtual

Check if an Particle has a collision in the grid; avoids multiple checks.

Test collisions between object and all objects in hgrid.

Implements HGRID_base.

Definition at line 83 of file HGRID_3D.cc.

References BOTTOMUP, HGrid::cellSizes_, CheckCell(), CheckCell_current(), BaseParticle::get_HGRID_Level(), BaseParticle::get_HGRID_x(), BaseParticle::get_HGRID_y(), BaseParticle::get_HGRID_z(), BaseParticle::get_InteractionRadius(), BaseParticle::get_Position(), HGRID_base::getHGridMethod(), HGrid::invCellSizes_, HGrid::occupiedLevelsMask, TOPDOWN, Vec3D::X, Vec3D::Y, and Vec3D::Z.

84 {
85  unsigned int startLevel = obj->get_HGRID_Level();
86 
87  int x, y, z;
88  Mdouble inv_size;
89 
90  switch(getHGridMethod())
91  {
92  case BOTTOMUP:
93  {
94  int occupiedLevelsMask = grid->occupiedLevelsMask >> obj->get_HGRID_Level();
95  for (unsigned int level = startLevel; level < grid->cellSizes_.size(); occupiedLevelsMask >>= 1, level++)
96  {
97  // If no objects in rest of grid, stop now
98  if (occupiedLevelsMask == 0) break;
99 
100  // If no objects at this level, go on to the next level
101  if ((occupiedLevelsMask & 1) == 0) continue;
102 
103  if (level == startLevel)
104  {
105  x = obj->get_HGRID_x();
106  y = obj->get_HGRID_y();
107  z = obj->get_HGRID_z();
108 
109  CheckCell_current(x, y, z, level, grid);
110 
111  CheckCell(x+1, y-1, z , level, obj, grid);
112  CheckCell(x+1, y, z , level, obj, grid);
113  CheckCell(x+1, y+1, z , level, obj, grid);
114 
115  CheckCell(x+1, y-1, z+1, level, obj, grid);
116  CheckCell(x+1, y, z+1, level, obj, grid);
117  CheckCell(x+1, y+1, z+1, level, obj, grid);
118 
119  CheckCell(x+1, y-1, z-1, level, obj, grid);
120  CheckCell(x+1, y, z-1, level, obj, grid);
121  CheckCell(x+1, y+1, z-1, level, obj, grid);
122 
123  CheckCell(x , y+1, z , level, obj, grid);
124  CheckCell(x , y, z-1, level, obj, grid);
125  CheckCell(x , y+1, z-1, level, obj, grid);
126  CheckCell(x , y+1, z+1, level, obj, grid);
127  }
128  else
129  {
130  int xs,ys,zs,xe,ye,ze;
131  inv_size = grid->invCellSizes_[level];
132  xs = (int) floor((obj->get_Position().X-obj->get_InteractionRadius()) * inv_size - 0.5);
133  xe = (int) floor((obj->get_Position().X+obj->get_InteractionRadius()) * inv_size + 0.5);
134  ys = (int) floor((obj->get_Position().Y-obj->get_InteractionRadius()) * inv_size - 0.5);
135  ye = (int) floor((obj->get_Position().Y+obj->get_InteractionRadius()) * inv_size + 0.5);
136  zs = (int) floor((obj->get_Position().Z-obj->get_InteractionRadius()) * inv_size - 0.5);
137  ze = (int) floor((obj->get_Position().Z+obj->get_InteractionRadius()) * inv_size + 0.5);
138  for(x=xs;x<=xe;x++)
139  {
140  for(y=ys;y<=ye;y++)
141  {
142  for(z=zs;z<=ze;z++)
143  {
144  CheckCell(x,y,z,level,obj,grid);
145  }
146  }
147  }
148  }
149  }
150  break;
151  }
152  case TOPDOWN:
153  {
154  int occupiedLevelsMask = grid->occupiedLevelsMask;
155  for (unsigned int level = 0; level <=startLevel; occupiedLevelsMask >>= 1, level++)
156  {
157  // If no objects in rest of grid, stop now
158  if (occupiedLevelsMask == 0) break;
159 
160  // If no objects at this level, go on to the next level
161  if ((occupiedLevelsMask & 1) == 0) continue;
162 
163  if (level == startLevel)
164  {
165  x = obj->get_HGRID_x();
166  y = obj->get_HGRID_y();
167  z = obj->get_HGRID_z();
168 
169  CheckCell_current(x, y, z, level, grid);
170 
171  CheckCell(x+1, y-1, z , level, obj, grid);
172  CheckCell(x+1, y, z , level, obj, grid);
173  CheckCell(x+1, y+1, z , level, obj, grid);
174 
175  CheckCell(x+1, y-1, z+1, level, obj, grid);
176  CheckCell(x+1, y, z+1, level, obj, grid);
177  CheckCell(x+1, y+1, z+1, level, obj, grid);
178 
179  CheckCell(x+1, y-1, z-1, level, obj, grid);
180  CheckCell(x+1, y, z-1, level, obj, grid);
181  CheckCell(x+1, y+1, z-1, level, obj, grid);
182 
183  CheckCell(x , y+1, z , level, obj, grid);
184  CheckCell(x , y, z-1, level, obj, grid);
185  CheckCell(x , y+1, z-1, level, obj, grid);
186  CheckCell(x , y+1, z+1, level, obj, grid);
187  }
188  else
189  {
190  int xs,ys,zs,xe,ye,ze;
191  inv_size = grid->invCellSizes_[level];
192  xs = (int) floor((obj->get_Position().X-obj->get_InteractionRadius()) * inv_size - 0.5);
193  xe = (int) floor((obj->get_Position().X+obj->get_InteractionRadius()) * inv_size + 0.5);
194  ys = (int) floor((obj->get_Position().Y-obj->get_InteractionRadius()) * inv_size - 0.5);
195  ye = (int) floor((obj->get_Position().Y+obj->get_InteractionRadius()) * inv_size + 0.5);
196  zs = (int) floor((obj->get_Position().Z-obj->get_InteractionRadius()) * inv_size - 0.5);
197  ze = (int) floor((obj->get_Position().Z+obj->get_InteractionRadius()) * inv_size + 0.5);
198  for(x=xs;x<=xe;x++)
199  {
200  for(y=ys;y<=ye;y++)
201  {
202  for(z=zs;z<=ze;z++)
203  {
204  CheckCell(x,y,z,level,obj,grid);
205  }
206  }
207  }
208  }
209  }
210  break;
211  }
212  }
213 }
HGridMethod getHGridMethod()
Definition: HGRID_base.cc:298
Mdouble X
Definition: Vector.h:44
Mdouble get_InteractionRadius() const
virtual void CheckCell(int x, int y, int z, int l, BaseParticle *obj, HGrid *grid)
Check collisions for a general cell.
Definition: HGRID_3D.cc:58
int get_HGRID_z() const
int get_HGRID_x() const
std::vector< double > invCellSizes_
Definition: HGRID.h:51
std::vector< double > cellSizes_
Definition: HGRID.h:50
double Mdouble
Definition: ExtendedMath.h:33
const Vec3D & get_Position() const
Mdouble Y
Definition: Vector.h:44
int get_HGRID_y() const
int occupiedLevelsMask
l-th bit of occupiedLevelsMask is 1 if level l is contains particles; initially zero (Implies max 32 ...
Definition: HGRID.h:55
virtual void CheckCell_current(int x, int y, int z, int l, HGrid *grid)
Checks for a collision in the particles own cell.
Definition: HGRID_3D.cc:31
Mdouble Z
Definition: Vector.h:44
int get_HGRID_Level() const
void HGRID_3D::CheckObjAgainstWholeGrid ( HGrid grid,
BaseParticle obj 
)
protected

Check if an Particle has a collision in the grid.

void HGRID_3D::constructor ( )
inline

This is the actually constructor it is called do both constructors above.

Definition at line 69 of file HGRID_3D.h.

References MD::set_dim(), and MD::set_dim_particle().

Referenced by HGRID_3D().

69  {
71  set_dim(3);
72  }
void set_dim_particle(int new_, unsigned int indSpecies=0)
Allows the dimension of the particle (f.e. for mass) to be changed.
Definition: MD.h:261
void set_dim(int new_dim)
Allows the dimension of the simulation to be changed.
Definition: MD.h:367
void HGRID_3D::HGRID_RemoveParticleFromHgrid ( BaseParticle obj)
protected

Definition at line 269 of file HGRID_3D.cc.

References HGrid::ComputeHashBucketIndex(), BaseParticle::get_HGRID_Level(), BaseParticle::get_HGRID_NextObject(), BaseParticle::get_HGRID_PrevObject(), BaseParticle::get_HGRID_x(), BaseParticle::get_HGRID_y(), BaseParticle::get_HGRID_z(), HGRID_base::grid, HGrid::objectBucket, BaseParticle::set_HGRID_NextObject(), and BaseParticle::set_HGRID_PrevObject().

Referenced by HGRID_UpdateParticleInHgrid().

270 {
271  int bucket = grid->ComputeHashBucketIndex(obj->get_HGRID_x(),obj->get_HGRID_y(), obj->get_HGRID_z(), obj->get_HGRID_Level());
272  if(obj->get_HGRID_PrevObject())
274  else
275  if(grid->objectBucket[bucket]==obj)
276  grid->objectBucket[bucket]=obj->get_HGRID_NextObject();
277 
278  if(obj->get_HGRID_NextObject())
280 }
void set_HGRID_PrevObject(BaseParticle *_new)
BaseParticle * get_HGRID_PrevObject() const
int get_HGRID_z() const
BaseParticle * get_HGRID_NextObject() const
int get_HGRID_x() const
int ComputeHashBucketIndex(int x, int y, int z, int l)
Computes hash bucket index in range [0, NUM_BUCKETS-1].
Definition: HGRID.cc:93
HGrid * grid
Definition: HGRID_base.h:131
int get_HGRID_y() const
void set_HGRID_NextObject(BaseParticle *_new)
int get_HGRID_Level() const
BaseParticle ** objectBucket
objectBucket[b] stores pointer to first element in hash bucket b; initially all NULL ...
Definition: HGRID.h:58
void HGRID_3D::HGRID_UpdateParticleInHgrid ( BaseParticle obj)
protected

This adds a partcile to the Grid, called in the grid setup routies.

Add object to grid square, and remember cell and level numbers, treating level as a third dimension coordinate.

Definition at line 220 of file HGRID_3D.cc.

References HGrid::ComputeHashBucketIndex(), BaseParticle::get_HGRID_Level(), BaseParticle::get_HGRID_x(), BaseParticle::get_HGRID_y(), BaseParticle::get_HGRID_z(), BaseParticle::get_PeriodicFromParticle(), BaseParticle::get_Position(), HGRID_base::grid, HGRID_RemoveParticleFromHgrid(), HGrid::invCellSizes_, HGrid::objectBucket, BaseParticle::set_HGRID_NextObject(), BaseParticle::set_HGRID_PrevObject(), BaseParticle::set_HGRID_x(), BaseParticle::set_HGRID_y(), BaseParticle::set_HGRID_z(), Vec3D::X, Vec3D::Y, and Vec3D::Z.

Referenced by Chute::add_particle(), and Chute::IsInsertable().

221 {
222  int l = obj->get_HGRID_Level();
223  Mdouble inv_size = grid->invCellSizes_[l];
224 
225  int x=(int)floor(obj->get_Position().X * inv_size);
226  int y=(int)floor(obj->get_Position().Y * inv_size);
227  int z=(int)floor(obj->get_Position().Z * inv_size);
228 
229 #ifdef ContactListHgrid
230  if(obj->get_HGRID_x()!=x||obj->get_HGRID_y()!=y||obj->get_HGRID_z()!=z||obj->get_PeriodicFromParticle())
231  {
232  if(!obj->get_PeriodicFromParticle())
233  {
234  int bucket = grid->ComputeHashBucketIndex(x, y, z, l);
235 
236  //First the object has to be removed
238  //Also remove all contact associated with it
239  getPossibleContactList().remove_ParticlePosibleContacts(obj);
240 
241  //And now reinserted
242  obj->set_HGRID_NextObject(grid->objectBucket[bucket]);
243  if(grid->objectBucket[bucket])
244  grid->objectBucket[bucket]->set_HGRID_PrevObject(obj);
245  obj->set_HGRID_PrevObject(0);
246  grid->objectBucket[bucket] = obj;
247  }
248 
249  obj->set_HGRID_x(x);
250  obj->set_HGRID_y(y);
251  obj->set_HGRID_z(z);
252  InsertObjAgainstGrid(grid, obj);
253  }
254 #else
255  int bucket = grid->ComputeHashBucketIndex(x,y,z,l);
256 
257  obj->set_HGRID_NextObject(grid->objectBucket[bucket]);
258  if(grid->objectBucket[bucket])
259  grid->objectBucket[bucket]->set_HGRID_PrevObject(obj);
260  obj->set_HGRID_PrevObject(0);
261  grid->objectBucket[bucket] = obj;
262 
263  obj->set_HGRID_x(x);
264  obj->set_HGRID_y(y);
265  obj->set_HGRID_z(z);
266 #endif
267 }
void set_HGRID_x(const int _new)
void set_HGRID_PrevObject(BaseParticle *_new)
Mdouble X
Definition: Vector.h:44
BaseParticle * get_PeriodicFromParticle() const
void set_HGRID_z(const int _new)
int get_HGRID_z() const
void HGRID_RemoveParticleFromHgrid(BaseParticle *obj)
Definition: HGRID_3D.cc:269
int get_HGRID_x() const
std::vector< double > invCellSizes_
Definition: HGRID.h:51
double Mdouble
Definition: ExtendedMath.h:33
void set_HGRID_y(const int _new)
int ComputeHashBucketIndex(int x, int y, int z, int l)
Computes hash bucket index in range [0, NUM_BUCKETS-1].
Definition: HGRID.cc:93
HGrid * grid
Definition: HGRID_base.h:131
const Vec3D & get_Position() const
Mdouble Y
Definition: Vector.h:44
int get_HGRID_y() const
void set_HGRID_NextObject(BaseParticle *_new)
Mdouble Z
Definition: Vector.h:44
int get_HGRID_Level() const
BaseParticle ** objectBucket
objectBucket[b] stores pointer to first element in hash bucket b; initially all NULL ...
Definition: HGRID.h:58
bool HGRID_3D::TestCell ( int  x,
int  y,
int  z,
int  l,
BaseParticle obj,
HGrid grid 
)
protected

Tests obj against all particles in cell similar to CheckCell, but links to TestObject instead of compute_internal_forces.

Definition at line 286 of file HGRID_3D.cc.

References HGrid::ComputeHashBucketIndex(), BaseParticle::get_HGRID_Level(), BaseParticle::get_HGRID_NextObject(), BaseParticle::get_HGRID_x(), BaseParticle::get_HGRID_y(), BaseParticle::get_HGRID_z(), HGrid::objectBucket, and HGRID_base::TestObject().

Referenced by TestObjAgainstGrid().

287 {
288  bool Test = true;
289 
290  // Loop through all objects in the bucket to find nearby objects
291  int bucket = grid->ComputeHashBucketIndex(x,y,z,l);
292  BaseParticle *p = grid->objectBucket[bucket];
293 
294  while (Test && p!=NULL)
295  {
296  if ((p->get_HGRID_x() == x) && (p->get_HGRID_y() == y) && (p->get_HGRID_z() == z) && (p->get_HGRID_Level() == l))
297  {
298  Test = Test && TestObject(obj,p);
299  }
300  p = p->get_HGRID_NextObject();
301  }
302 
303  return Test;
304 }
int get_HGRID_z() const
BaseParticle * get_HGRID_NextObject() const
int get_HGRID_x() const
virtual bool TestObject(BaseParticle *pI, BaseParticle *pJ)
criterium for inserting a particle (returns false, if particles overlap;)
Definition: HGRID_base.cc:251
int ComputeHashBucketIndex(int x, int y, int z, int l)
Computes hash bucket index in range [0, NUM_BUCKETS-1].
Definition: HGRID.cc:93
int get_HGRID_y() const
int get_HGRID_Level() const
BaseParticle ** objectBucket
objectBucket[b] stores pointer to first element in hash bucket b; initially all NULL ...
Definition: HGRID.h:58
bool HGRID_3D::TestObjAgainstGrid ( HGrid grid,
BaseParticle obj 
)
protected

Tests obj against all neighbouring particles similar to CheckObjAgainstGrid, but links to TestCell instead of CheckCell.

Definition at line 310 of file HGRID_3D.cc.

References HGrid::cellSizes_, BaseParticle::get_Position(), HGrid::invCellSizes_, HGrid::occupiedLevelsMask, TestCell(), Vec3D::X, Vec3D::Y, and Vec3D::Z.

Referenced by Chute::IsInsertable().

311 {
312  bool Test = true;
313 
314  int x, y, z;
315  Mdouble inv_size;
316  int occupiedLevelsMask = grid->occupiedLevelsMask;
317 
318  for (unsigned int level = 0; level < grid->cellSizes_.size() && Test; occupiedLevelsMask >>= 1, level++)
319  {
320  // If no objects in rest of grid, stop now
321  if (occupiedLevelsMask == 0) break;
322 
323  // If no objects at this level, go on to the next level
324  if ((occupiedLevelsMask & 1) == 0) continue;
325 
326  // Treat level as a third dimension coordinate
327  inv_size = grid->invCellSizes_[level];
328 
329  x = (int)floor(obj->get_Position().X * inv_size);
330  y = (int)floor(obj->get_Position().Y * inv_size);
331  z = (int)floor(obj->get_Position().Z * inv_size);
332 
333  Test = Test
334  && TestCell(x , y , z , level, obj, grid)
335  && TestCell(x , y-1, z , level, obj, grid)
336  && TestCell(x , y+1, z , level, obj, grid)
337  && TestCell(x-1, y , z , level, obj, grid)
338  && TestCell(x+1, y , z , level, obj, grid)
339  && TestCell(x-1, y-1, z , level, obj, grid)
340  && TestCell(x-1, y+1, z , level, obj, grid)
341  && TestCell(x+1, y-1, z , level, obj, grid)
342  && TestCell(x+1, y+1, z , level, obj, grid)
343 
344  && TestCell(x , y , z-1, level, obj, grid)
345  && TestCell(x , y-1, z-1, level, obj, grid)
346  && TestCell(x , y+1, z-1, level, obj, grid)
347  && TestCell(x-1, y , z-1, level, obj, grid)
348  && TestCell(x+1, y , z-1, level, obj, grid)
349  && TestCell(x-1, y-1, z-1, level, obj, grid)
350  && TestCell(x-1, y+1, z-1, level, obj, grid)
351  && TestCell(x+1, y-1, z-1, level, obj, grid)
352  && TestCell(x+1, y+1, z-1, level, obj, grid)
353 
354  && TestCell(x , y , z+1, level, obj, grid)
355  && TestCell(x , y-1, z+1, level, obj, grid)
356  && TestCell(x , y+1, z+1, level, obj, grid)
357  && TestCell(x-1, y , z+1, level, obj, grid)
358  && TestCell(x+1, y , z+1, level, obj, grid)
359  && TestCell(x-1, y-1, z+1, level, obj, grid)
360  && TestCell(x-1, y+1, z+1, level, obj, grid)
361  && TestCell(x+1, y-1, z+1, level, obj, grid)
362  && TestCell(x+1, y+1, z+1, level, obj, grid);
363  } //end for level
364 
365  return Test;
366 }
Mdouble X
Definition: Vector.h:44
std::vector< double > invCellSizes_
Definition: HGRID.h:51
std::vector< double > cellSizes_
Definition: HGRID.h:50
bool TestCell(int x, int y, int z, int l, BaseParticle *obj, HGrid *grid)
Tests obj against all particles in cell similar to CheckCell, but links to TestObject instead of comp...
Definition: HGRID_3D.cc:286
double Mdouble
Definition: ExtendedMath.h:33
const Vec3D & get_Position() const
Mdouble Y
Definition: Vector.h:44
int occupiedLevelsMask
l-th bit of occupiedLevelsMask is 1 if level l is contains particles; initially zero (Implies max 32 ...
Definition: HGRID.h:55
Mdouble Z
Definition: Vector.h:44

The documentation for this class was generated from the following files: