Description
Classes | |
struct | chrono::granular::ChGranSphereData |
Holds pointers to kinematic quantities of the granular system. These pointers must be in device-accessible memory. More... | |
struct | chrono::granular::ChGranParams |
Parameters needed for sphere-based granular dynamics. More... | |
class | chrono::granular::ChSystemGranularSMC |
Main Chrono::Granular system class used to control and dispatch the GPU sphere-only solver. More... | |
struct | chrono::granular::ChTriangleSoup< T3 > |
Class used to hold pointers for mesh arrays. More... | |
struct | chrono::granular::ChGranMeshFamilyFrame< T > |
Position and rotation matrix defining the frame of a triangle mesh. More... | |
struct | chrono::granular::ChGranParams_trimesh |
Extra parameters needed for triangle-sphere contact. More... | |
class | chrono::granular::ChSystemGranularSMC_trimesh |
Class implements functionality required to handle the interaction between a mesh soup and granular material. More... | |
Macros | |
#define | GET_OUTPUT_SETTING(setting) (this->output_flags & setting) |
Typedefs | |
typedef std::function< double3(float)> | chrono::granular::GranPositionFunction |
Used to compute position as a function of time. | |
Enumerations | |
enum | chrono::granular::GRAN_VERBOSITY { QUIET = 0, INFO = 1, METRICS = 2 } |
Verbosity level of the system. | |
enum | chrono::granular::GRAN_OUTPUT_MODE { CSV, BINARY, HDF5, NONE } |
Output mode of system. | |
enum | chrono::granular::GRAN_TIME_INTEGRATOR { FORWARD_EULER, CHUNG, CENTERED_DIFFERENCE, EXTENDED_TAYLOR } |
How are we integrating through time. | |
enum | chrono::granular::GRAN_FRICTION_MODE { FRICTIONLESS, SINGLE_STEP, MULTI_STEP } |
Supported friction model. | |
enum | chrono::granular::GRAN_ROLLING_MODE { NO_RESISTANCE, SCHWARTZ, ELASTIC_PLASTIC } |
Rolling resistance models – ELASTIC_PLASTIC not implemented yet. | |
enum | GRAN_OUTPUT_FLAGS { ABSV = 1, VEL_COMPONENTS = 2, FIXITY = 4, ANG_VEL_COMPONENTS = 8 } |
Functions | |
__device__ __host__ int64_t3 | convertPosLocalToGlobal (unsigned int ownerSD, const int3 &local_pos, GranParamsPtr gran_params) |
Convert position from its owner subdomain local frame to the global big domain frame. | |
__device__ void | figureOutTouchedSD (int64_t sphCenter_X_relative, int64_t sphCenter_Y_relative, int64_t sphCenter_Z_relative, unsigned int SDs[MAX_SDs_TOUCHED_BY_SPHERE], GranParamsPtr gran_params) |
Takes in a sphere's position and inserts into the given int array[8] which subdomains, if any, are touched The array is indexed with the ones bit equal to +/- x, twos bit equal to +/- y, and the fours bit equal to +/- z A bit set to 0 means the lower index, whereas 1 means the higher index (lower + 1) The kernel computes global x, y, and z indices for the bottom-left subdomain and then uses those to figure out which subdomains described in the corresponding 8-SD cube are touched by the sphere. More... | |
template<unsigned int CUB_THREADS> | |
__global__ void | sphereBroadphase_dryrun (GranSphereDataPtr sphere_data, unsigned int nSpheres, GranParamsPtr gran_params) |
This kernel call prepares information that will be used in a subsequent kernel that performs the actual time stepping. More... | |
template<unsigned int CUB_THREADS> | |
__global__ void | sphereBroadphase (GranSphereDataPtr sphere_data, unsigned int nSpheres, GranParamsPtr gran_params) |
__device__ int3 | getOffsetFromSDs (unsigned int thisSD, unsigned int otherSD, GranParamsPtr gran_params) |
Get position offset between two SDs. | |
__device__ void | findNewLocalCoords (GranSphereDataPtr sphere_data, unsigned int mySphereID, int64_t global_pos_X, int64_t global_pos_Y, int64_t global_pos_Z, GranParamsPtr gran_params) |
update local positions and SD based on global position | |
__device__ void | applyGravity (float3 &sphere_force, GranParamsPtr gran_params) |
__device__ void | applyExternalForces_frictionless (unsigned int ownerSD, const int3 &sphPos_local, const float3 &sphVel, float3 &sphere_force, GranParamsPtr gran_params, GranSphereDataPtr sphere_data, BC_type *bc_type_list, BC_params_t< int64_t, int64_t3 > *bc_params_list, unsigned int nBCs) |
Compute forces on a sphere from walls, BCs, and gravity. | |
__device__ void | applyExternalForces (unsigned int currSphereID, unsigned int ownerSD, const int3 &sphPos_local, const float3 &sphVel, const float3 &sphOmega, float3 &sphere_force, float3 &sphere_ang_acc, GranParamsPtr gran_params, GranSphereDataPtr sphere_data, BC_type *bc_type_list, BC_params_t< int64_t, int64_t3 > *bc_params_list, unsigned int nBCs) |
Compute forces on a sphere from walls, BCs, and gravity. | |
__device__ float3 | computeSphereNormalForces (float &reciplength, float3 &vrel_t, float3 &delta_r, const int3 &sphereA_pos, const int3 &sphereB_pos, const float3 &sphereA_vel, const float3 &sphereB_vel, GranParamsPtr gran_params) |
Compute normal forces for a contacting pair. | |
__device__ float | integrateForwardEuler (float stepsize_SU, float val_dt) |
Compute update for a quantity using Forward Euler integrator. | |
__device__ float | integrateChung_vel (float stepsize_SU, float acc, float acc_old) |
Compute update for a velocity using Chung integrator. | |
__device__ float | integrateChung_pos (float stepsize_SU, float vel_old, float acc, float acc_old) |
Compute update for a position using Chung integrator. | |
Variables | |
const GranPositionFunction | chrono::granular::GranPosFunction_default = [](float t) { return make_double3(0, 0, 0); } |
Position function representing no motion or offset as a funtion of time. | |
Function Documentation
◆ figureOutTouchedSD()
|
inline |
Takes in a sphere's position and inserts into the given int array[8] which subdomains, if any, are touched The array is indexed with the ones bit equal to +/- x, twos bit equal to +/- y, and the fours bit equal to +/- z A bit set to 0 means the lower index, whereas 1 means the higher index (lower + 1) The kernel computes global x, y, and z indices for the bottom-left subdomain and then uses those to figure out which subdomains described in the corresponding 8-SD cube are touched by the sphere.
The kernel then converts these indices to indices into the global SD list via the (currently local) conv[3] data structure Should be mostly bug-free, especially away from boundaries
◆ sphereBroadphase()
__global__ void sphereBroadphase | ( | GranSphereDataPtr | sphere_data, |
unsigned int | nSpheres, | ||
GranParamsPtr | gran_params | ||
) |
Set aside shared memory
◆ sphereBroadphase_dryrun()
__global__ void sphereBroadphase_dryrun | ( | GranSphereDataPtr | sphere_data, |
unsigned int | nSpheres, | ||
GranParamsPtr | gran_params | ||
) |
This kernel call prepares information that will be used in a subsequent kernel that performs the actual time stepping.
Template arguments:
- CUB_THREADS: the number of threads used in this kernel, comes into play when invoking CUB block collectives
Assumptions:
- Granular material is made up of monodisperse spheres.
- The function below assumes the spheres are in a box
- The box has dimensions L x D x H.
- The reference frame associated with the box:
- The x-axis is along the length L of the box
- The y-axis is along the width D of the box
- The z-axis is along the height H of the box
- A sphere cannot touch more than eight SDs
Basic idea: use domain decomposition on the rectangular box and figure out how many SDs each sphere touches. The subdomains are axis-aligned relative to the reference frame associated with the box. The origin of the box is at the center of the box. The orientation of the box is defined relative to a world inertial reference frame.
Nomenclature:
- SD: subdomain.
- BD: the big-domain, which is the union of all SDs
- NULL_GRANULAR_ID: the equivalent of a non-sphere SD ID, or a non-sphere ID
Notes:
- The SD with ID=0 is the catch-all SD. This is the SD in which a sphere ends up if its not inside the rectangular box. Usually, there is no sphere in this SD (THIS IS NOT IMPLEMENTED AS SUCH FOR NOW)
Set aside shared memory