Description

Main Chrono::Granular system class used to control and dispatch the GPU sphere-only solver.

#include <ChGranular.h>

Inheritance diagram for chrono::granular::ChSystemGranularSMC:
Collaboration diagram for chrono::granular::ChSystemGranularSMC:

Public Member Functions

 ChSystemGranularSMC (float sphere_rad, float density, float3 boxDims)
 Construct granular system with given sphere radius, density, and big domain dimensions.
 
unsigned int get_SD_count () const
 Return number of subdomains in the big domain.
 
void set_gravitational_acceleration (float xVal, float yVal, float zVal)
 Set acceleration of gravity.
 
size_t getNumSpheres () const
 Return total number of spheres in the system.
 
size_t estimateMemUsage () const
 Roughly estimates the total amount of memory used by the system.
 
size_t Create_BC_Sphere (float center[3], float radius, bool outward_normal, bool track_forces)
 Create an axis-aligned sphere boundary condition.
 
size_t Create_BC_Cone_Z (float cone_tip[3], float slope, float hmax, float hmin, bool outward_normal, bool track_forces)
 Create an z-axis aligned cone boundary condition.
 
size_t Create_BC_Plane (float plane_pos[3], float plane_normal[3], bool track_forces)
 Create plane boundary condition.
 
size_t Create_BC_Cyl_Z (float center[3], float radius, bool outward_normal, bool track_forces)
 Create an z-axis aligned cylinder boundary condition.
 
void createWallBCs ()
 Create big domain walls out of plane boundary conditions.
 
bool disable_BC_by_ID (size_t BC_id)
 Disable a boundary condition by its ID, returns false if the BC does not exist.
 
bool enable_BC_by_ID (size_t BC_id)
 Enable a boundary condition by its ID, returns false if the BC does not exist.
 
bool set_BC_offset_function (size_t BC_id, const GranPositionFunction &offset_function)
 Enable a boundary condition by its ID, returns false if the BC does not exist.
 
bool getBCReactionForces (size_t BC_id, float forces[3]) const
 Get the reaction forces on a boundary by ID, returns false if the forces are invalid (bad BC ID)
 
void setOutputMode (GRAN_OUTPUT_MODE mode)
 Set the output mode of the simulation.
 
void setVerbose (GRAN_VERBOSITY level)
 Set simualtion verbosity – used to check on very large, slow simulations or debug.
 
void setOutputFlags (unsigned char flags)
 Set output settings bit flags by bitwise ORing settings in GRAN_OUTPUT_FLAGS.
 
void set_fixed_stepSize (float size_UU)
 Set timestep size.
 
void disableMinLength ()
 Ensure that the deformation-based length unit is used.
 
void set_timeIntegrator (GRAN_TIME_INTEGRATOR new_integrator)
 Set the time integration scheme for the system.
 
void set_friction_mode (GRAN_FRICTION_MODE new_mode)
 Set friction formulation. More...
 
void set_rolling_mode (GRAN_ROLLING_MODE new_mode)
 Set rolling resistence formulation. More...
 
double get_max_z () const
 Get the max z position of the spheres, allows easier co-simulation.
 
virtual double advance_simulation (float duration)
 Advance simulation by duration in user units, return actual duration elapsed Requires initialize() to have been called.
 
virtual void initialize ()
 Initialize simulation so that it can be advanced Must be called before advance_simulation and after simulation parameters are set.
 
void set_static_friction_coeff_SPH2SPH (float mu)
 Set sphere-to-sphere static friction coefficient.
 
void set_static_friction_coeff_SPH2WALL (float mu)
 Set sphere-to-wall static friction coefficient.
 
void set_rolling_coeff_SPH2SPH (float mu)
 Set sphere-to-sphere rolling friction coefficient – units and use vary by rolling friction mode.
 
void set_rolling_coeff_SPH2WALL (float mu)
 Set sphere-to-wall rolling friction coefficient – units and use vary by rolling friction mode.
 
void set_spinning_coeff_SPH2SPH (float mu)
 Set sphere-to-sphere spinning friction coefficient – units and use vary by spinning friction mode.
 
void set_spinning_coeff_SPH2WALL (float mu)
 Set sphere-to-wall spinning friction coefficient – units and use vary by spinning friction mode.
 
void set_K_n_SPH2SPH (double someValue)
 Set sphere-to-sphere normal contact stiffness.
 
void set_K_n_SPH2WALL (double someValue)
 Set sphere-to-wall normal contact stiffness.
 
void set_Gamma_n_SPH2SPH (double someValue)
 Set sphere-to-sphere normal damping coefficient.
 
void set_Gamma_n_SPH2WALL (double someValue)
 Set sphere-to-wall normal damping coefficient.
 
void set_K_t_SPH2SPH (double someValue)
 Set sphere-to-sphere tangent contact stiffness.
 
void set_Gamma_t_SPH2SPH (double someValue)
 Set sphere-to-sphere tangent damping coefficient.
 
void set_K_t_SPH2WALL (double someValue)
 Set sphere-to-wall tangent contact stiffness.
 
void set_Gamma_t_SPH2WALL (double someValue)
 Set sphere-to-wall tangent damping coefficient.
 
void set_Cohesion_ratio (float someValue)
 Set the ratio of cohesion to gravity for monodisperse spheres. Assumes a constant cohesion model.
 
void set_Adhesion_ratio_S2W (float someValue)
 Set the ratio of adhesion to gravity for sphere to wall. Assumes a constant cohesion model.
 
void set_BD_Fixed (bool fixed)
 Set the big domain to be fixed or not; if fixed it will ignore any given position functions.
 
void setParticlePositions (const std::vector< float3 > &points, const std::vector< float3 > &vels=std::vector< float3 >(), const std::vector< float3 > &ang_vels=std::vector< float3 >())
 Set initial particle positions. MUST be called only once and MUST be called before initialize.
 
void setParticleFixed (const std::vector< bool > &fixed)
 Set particle fixity. MUST be called only once and MUST be called before initialize.
 
float3 getPosition (int nSphere)
 return particle position given sphere index
 
float getAbsVelocity (int nSphere)
 
float3 getVelocity (int nSphere)
 
float3 getAngularVelocity (int nSphere)
 
int getNumContacts ()
 
float3 Get_BC_Plane_Position (size_t plane_id)
 
void setBDWallsMotionFunction (const GranPositionFunction &pos_fn)
 Prescribe the motion of the big domain, allows wavetank-style simulations.
 
void setRecordingContactInfo (bool record)
 
void setPsiFactors (unsigned int psi_T_new, unsigned int psi_L_new, float psi_R_new=1.f)
 Set tuning psi factors for tuning the non-dimensionalization.
 
void checkSDCounts (std::string ofile, bool write_out, bool verbose) const
 Copy back the subdomain device data and save it to a file for error checking on the priming kernel.
 
void writeFile (std::string ofile) const
 Writes out particle positions according to the system output mode.
 
void writeContactInfoFile (std::string ofile) const
 Writes out contact info.
 
void setMaxSafeVelocity_SU (float max_vel)
 Safety check velocity to ensure the simulation is still stable.
 

Public Attributes

GranPositionFunction BDOffsetFunction
 The offset function for the big domain walls.
 

Protected Member Functions

int3 getSDTripletFromID (unsigned int SD_ID) const
 Wrap the device helper function.
 
void initializeSpheres ()
 Create a helper to do sphere initialization.
 
void defragment_initial_positions ()
 Sorts particle positions spatially in order to improve memory efficiency.
 
void setupSphereDataStructures ()
 Setup sphere data, initialize local coords.
 
void partitionBD ()
 This method figures out how big a SD is, and how many SDs are going to be necessary in order to cover the entire BD. More...
 
void copyConstSphereDataToDevice ()
 Copy constant sphere data to device.
 
void resetBroadphaseInformation ()
 Reset binning and broadphase info.
 
void resetSphereAccelerations ()
 Reset sphere accelerations.
 
void resetBCForces ()
 Reset sphere-wall forces.
 
void packSphereDataPointers ()
 Collect all the sphere data into the member struct.
 
void runSphereBroadphase ()
 Run the first sphere broadphase pass to get things started.
 
template<typename T1 , typename T2 >
T1 convertToPosSU (T2 val)
 Helper function to convert a position in UU to its SU representation while also changing data type.
 
void convertBCUnits ()
 convert all BCs from UU to SU
 
float get_max_vel () const
 Max velocity of all particles in system.
 
virtual double get_max_K () const
 Get the maximum stiffness term in the system.
 
virtual void switchToSimUnits ()
 This method defines the mass, time, length Simulation Units. More...
 
void setBCOffset (const BC_type &, const BC_params_t< float, float3 > &params_UU, BC_params_t< int64_t, int64_t3 > &params_SU, double3 offset_UU)
 Set the position function of a boundary condition and account for the offset.
 
void updateBCPositions ()
 Update positions of each boundary condition using prescribed functions.
 

Protected Attributes

double LENGTH_SU2UU
 1 / C_L. Any length expressed in SU is a multiple of SU2UU
 
double TIME_SU2UU
 1 / C_T. Any time quantity in SU is measured as a positive multiple of TIME_SU2UU
 
double MASS_SU2UU
 1 / C_M. Any mass quantity is measured as a positive multiple of MASS_SU2UU.
 
double FORCE_SU2UU
 1 / C_F. Any force quantity is measured as a multiple of FORCE_SU2UU.
 
double TORQUE_SU2UU
 1 / C_tau. Any torque quantity is measured as a multiple of TORQUE_SU2UU.
 
double VEL_SU2UU
 1 / C_v. Any velocity quantity is measured as a multiple of VEL_SU2UU.
 
unsigned int psi_T
 Safety factor on simulation time.
 
unsigned int psi_L
 Safety factor on space adim.
 
float psi_R
 Fraction of sphere radius which gives an upper bound on the length unit.
 
ChGranParamsgran_params
 Holds the sphere and big-domain-related params in unified memory.
 
ChGranSphereDatasphere_data
 Holds system degrees of freedom.
 
GRAN_VERBOSITY verbosity
 Allows the code to be very verbose for debugging.
 
bool use_min_length_unit
 If dividing the longest box dimension into INT_MAX pieces gives better resolution than the deformation-based scaling, do that.
 
unsigned char output_flags
 Bit flags indicating what fields to write out during writeFile Set with the GRAN_OUTPUT_FLAGS enum.
 
GRAN_OUTPUT_MODE file_write_mode
 How to write the output files? Default is CSV.
 
unsigned int nSpheres
 Number of discrete elements.
 
unsigned int nSDs
 Number of subdomains.
 
std::vector< int, cudallocator< int > > sphere_local_pos_X
 Store X positions relative to owner subdomain in unified memory.
 
std::vector< int, cudallocator< int > > sphere_local_pos_Y
 Store Y positions relative to owner subdomain in unified memory.
 
std::vector< int, cudallocator< int > > sphere_local_pos_Z
 Store Z positions relative to owner subdomain in unified memory.
 
std::vector< float, cudallocator< float > > pos_X_dt
 Store X velocity in unified memory.
 
std::vector< float, cudallocator< float > > pos_Y_dt
 Store Y velocity in unified memory.
 
std::vector< float, cudallocator< float > > pos_Z_dt
 Store Z velocity in unified memory.
 
std::vector< float, cudallocator< float > > sphere_Omega_X
 Store X angular velocity (axis and magnitude in one vector) in unified memory.
 
std::vector< float, cudallocator< float > > sphere_Omega_Y
 Store Y angular velocity (axis and magnitude in one vector) in unified memory.
 
std::vector< float, cudallocator< float > > sphere_Omega_Z
 Store Z angular velocity (axis and magnitude in one vector) in unified memory.
 
std::vector< float, cudallocator< float > > sphere_acc_X
 Store X acceleration in unified memory.
 
std::vector< float, cudallocator< float > > sphere_acc_Y
 Store Y acceleration in unified memory.
 
std::vector< float, cudallocator< float > > sphere_acc_Z
 Store Z acceleration in unified memory.
 
std::vector< float, cudallocator< float > > sphere_ang_acc_X
 Store X angular acceleration (axis and magnitude in one vector) in unified memory.
 
std::vector< float, cudallocator< float > > sphere_ang_acc_Y
 Store Y angular acceleration (axis and magnitude in one vector) in unified memory.
 
std::vector< float, cudallocator< float > > sphere_ang_acc_Z
 Store Z angular acceleration (axis and magnitude in one vector) in unified memory.
 
std::vector< float, cudallocator< float > > sphere_acc_X_old
 X acceleration history used for the Chung integrator in unified memory.
 
std::vector< float, cudallocator< float > > sphere_acc_Y_old
 Y acceleration history used for the Chung integrator in unified memory.
 
std::vector< float, cudallocator< float > > sphere_acc_Z_old
 Z acceleration history used for the Chung integrator in unified memory.
 
std::vector< float, cudallocator< float > > sphere_ang_acc_X_old
 X angular acceleration history used for the Chung integrator in unified memory.
 
std::vector< float, cudallocator< float > > sphere_ang_acc_Y_old
 Y angular acceleration history used for the Chung integrator in unified memory.
 
std::vector< float, cudallocator< float > > sphere_ang_acc_Z_old
 Z angular acceleration history used for the Chung integrator in unified memory.
 
std::vector< not_stupid_bool, cudallocator< not_stupid_bool > > sphere_fixed
 Fixity of each sphere.
 
std::vector< unsigned int, cudallocator< unsigned int > > contact_partners_map
 Set of contact partners for each sphere. Only used in frictional simulations.
 
std::vector< not_stupid_bool, cudallocator< not_stupid_bool > > contact_active_map
 Whether the frictional contact at an index is active.
 
std::vector< float3, cudallocator< float3 > > contact_history_map
 Tracks the tangential history vector for a given contact pair. Only used in multistep friction.
 
std::vector< float3, cudallocator< float3 > > normal_contact_force
 Tracks the normal contact force for a given contact pair.
 
std::vector< float3, cudallocator< float3 > > tangential_friction_force
 Tracks the tangential contact force for a given contact pair.
 
std::vector< float3, cudallocator< float3 > > rolling_friction_torque
 Tracks the rolling resistance for a given contact pair.
 
float X_accGrav
 X gravity in user units.
 
float Y_accGrav
 Y gravity in user units.
 
float Z_accGrav
 Z gravity in user units.
 
std::vector< unsigned int, cudallocator< unsigned int > > SD_NumSpheresTouching
 Entry "i" says how many spheres touch subdomain i.
 
std::vector< unsigned int, cudallocator< unsigned int > > SD_SphereCompositeOffsets
 
std::vector< unsigned int, cudallocator< unsigned int > > spheres_in_SD_composite
 Array containing the IDs of the spheres stored in the SDs associated with the box.
 
std::vector< unsigned int, cudallocator< unsigned int > > sphere_owner_SDs
 List of owner subdomains for each sphere.
 
float stepSize_UU
 User provided timestep in UU.
 
float stepSize_SU
 SU (adimensional) value of time step.
 
GRAN_TIME_INTEGRATOR time_integrator
 how is the system integrated through time?
 
float elapsedSimTime
 Total time elapsed since beginning of simulation.
 
double K_n_s2s_UU
 Normal stiffness for sphere-to-sphere contact (Hertzian spring constant)
 
double K_n_s2w_UU
 Normal stiffness for sphere-to-wall contact (Hertzian spring constant)
 
double Gamma_n_s2s_UU
 Normal damping for sphere-to-sphere.
 
double Gamma_n_s2w_UU
 Normal damping for sphere-to-wall.
 
double K_t_s2s_UU
 Tangential stiffness for sphere-to-sphere (Hertzian spring constant.
 
double K_t_s2w_UU
 Tangential stiffness for sphere-to-wall (Hertzian spring constant.
 
double Gamma_t_s2s_UU
 Tangential damping for sphere-to-sphere.
 
double Gamma_t_s2w_UU
 Tangential damping for sphere-to-wall.
 
double rolling_coeff_s2s_UU
 Rolling friction coefficient for sphere-to-sphere.
 
double rolling_coeff_s2w_UU
 Rolling friction coefficient for sphere-to-wall Units and use are dependent on the rolling friction model used.
 
double spinning_coeff_s2s_UU
 Spinning friction coefficient for sphere-to-sphere.
 
double spinning_coeff_s2w_UU
 Spinning friction coefficient for sphere-to-wall Units and use are dependent on the spinning friction model used.
 
float cohesion_over_gravity
 Store the ratio of the acceleration due to cohesion vs the acceleration due to gravity, makes simple API.
 
float adhesion_s2w_over_gravity
 Store the ratio of the acceleration due to adhesion vs the acceleration due to gravity.
 
int64_t3 BD_rest_frame_SU
 The reference point for UU to SU local coordinates.
 
std::vector< BC_type, cudallocator< BC_type > > BC_type_list
 List of generalized boundary conditions that constrain sphere motion.
 
std::vector< BC_params_t< int64_t, int64_t3 >, cudallocator< BC_params_t< int64_t, int64_t3 > > > BC_params_list_SU
 Sim-unit (adimensional) details of boundary conditions.
 
std::vector< BC_params_t< float, float3 >, cudallocator< BC_params_t< float, float3 > > > BC_params_list_UU
 User-unit (dimensional) details of boundary conditions.
 
std::vector< GranPositionFunctionBC_offset_function_list
 Offset motions functions for boundary conditions – used for moving walls, wavetanks, etc.
 
const float sphere_radius_UU
 User defined radius of the sphere.
 
const float sphere_density_UU
 User defined density of the sphere.
 
const float box_size_X
 X-length of the big domain; defines the global X axis located at the CM of the box.
 
const float box_size_Y
 Y-length of the big domain; defines the global Y axis located at the CM of the box.
 
const float box_size_Z
 Z-length of the big domain; defines the global Z axis located at the CM of the box.
 
std::vector< float3 > user_sphere_positions
 User-provided sphere positions in UU.
 
std::vector< float3 > user_sphere_vel
 User-provided sphere velocities in UU.
 
std::vector< float3 > user_sphere_ang_vel
 
std::vector< bool > user_sphere_fixed
 User-provided sphere fixity as bools.
 
bool BD_is_fixed = true
 Allow the user to set the big domain to be fixed, ignoring any given position functions.
 

Member Function Documentation

◆ partitionBD()

void chrono::granular::ChSystemGranularSMC::partitionBD ( )
protected

This method figures out how big a SD is, and how many SDs are going to be necessary in order to cover the entire BD.

Nomenclature: BD: Big domain, SD: Sub-domain.

◆ set_friction_mode()

void chrono::granular::ChSystemGranularSMC::set_friction_mode ( GRAN_FRICTION_MODE  new_mode)
inline

Set friction formulation.

The frictionless setting uses a streamlined solver and avoids storing any physics information associated with friction

◆ set_rolling_mode()

void chrono::granular::ChSystemGranularSMC::set_rolling_mode ( GRAN_ROLLING_MODE  new_mode)
inline

Set rolling resistence formulation.

NOTE: This requires friction to be active, otherwise this setting will be ignored

◆ switchToSimUnits()

void chrono::granular::ChSystemGranularSMC::switchToSimUnits ( )
protectedvirtual

This method defines the mass, time, length Simulation Units.

It also sets several other constants that enter the scaling of various physical quantities set by the user.

SU values for normal stiffnesses for S2S and S2W


The documentation for this class was generated from the following files:
  • /builds/uwsbel/chrono/src/chrono_granular/physics/ChGranular.h
  • /builds/uwsbel/chrono/src/chrono_granular/physics/ChGranular.cpp