Description
Functions | |
__device__ __host__ int64_t3 | convertPosLocalToGlobal (unsigned int ownerSD, const int3 &local_pos, ChSystemGpu_impl::GranParamsPtr gran_params) |
Convert position from its owner subdomain local frame to the global big domain frame. | |
__device__ void | figureOutTouchedSD (int sphCenter_X_local, int sphCenter_Y_local, int sphCenter_Z_local, int3 ownerSD, unsigned int SDs[MAX_SDs_TOUCHED_BY_SPHERE], ChSystemGpu_impl::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 | getNumberOfSpheresTouchingEachSD (ChSystemGpu_impl::GranSphereDataPtr sphere_data, unsigned int nSpheres, ChSystemGpu_impl::GranParamsPtr gran_params) |
Template arguments: More... | |
__device__ int3 | getOffsetFromSDs (unsigned int thisSD, unsigned int otherSD, ChSystemGpu_impl::GranParamsPtr gran_params) |
Get position offset between two SDs. | |
__device__ void | findNewLocalCoords (ChSystemGpu_impl::GranSphereDataPtr sphere_data, unsigned int mySphereID, int64_t global_pos_X, int64_t global_pos_Y, int64_t global_pos_Z, ChSystemGpu_impl::GranParamsPtr gran_params) |
update local positions and SD based on global position | |
__device__ void | applyGravity (float3 &sphere_force, ChSystemGpu_impl::GranParamsPtr gran_params) |
Apply gravity to a sphere. | |
__device__ void | applyExternalForces_frictionless (unsigned int ownerSD, const int3 &sphPos_local, const float3 &sphVel, float3 &sphere_force, ChSystemGpu_impl::GranParamsPtr gran_params, ChSystemGpu_impl::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, ChSystemGpu_impl::GranParamsPtr gran_params, ChSystemGpu_impl::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, ChSystemGpu_impl::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. | |
template<class IN_T , class IN_T3 , class OUT_T3 = IN_T3> | |
__device__ OUT_T3 | apply_frame_transform (const IN_T3 &point, const IN_T *pos, const IN_T *rot_mat) |
Point is in the LRF, rot_mat rotates LRF to GRF, pos translates LRF to GRF LRF: local reference frame GRF: global reference frame. | |
template<class T3 > | |
__device__ void | convert_pos_UU2SU (T3 &pos, ChSystemGpu_impl::GranParamsPtr gran_params) |
Convert position vector from user units to scaled units. | |
__inline__ __device__ void | triangle_figureOutSDBox (const float3 &vA, const float3 &vB, const float3 &vC, int *L, int *U, ChSystemGpu_impl::GranParamsPtr gran_params) |
Takes in a triangle ID and figures out an SD AABB for broadphase use. | |
__device__ unsigned int | triangle_countTouchedSDs (unsigned int triangleID, const ChSystemGpuMesh_impl::TriangleSoupPtr triangleSoup, ChSystemGpu_impl::GranParamsPtr gran_params, ChSystemGpuMesh_impl::MeshParamsPtr tri_params) |
Takes in a triangle's position in UU and finds out how many SDs it touches. More... | |
__device__ void | triangle_figureOutTouchedSDs (unsigned int triangleID, const ChSystemGpuMesh_impl::TriangleSoupPtr triangleSoup, unsigned int *touchedSDs, const ChSystemGpu_impl::GranParamsPtr &gran_params, const ChSystemGpuMesh_impl::MeshParamsPtr &tri_params) |
Takes in a triangle's position in UU and finds out what SDs it touches. More... | |
__global__ void | determineCountOfSDsTouchedByEachTriangle (const ChSystemGpuMesh_impl::TriangleSoupPtr d_triangleSoup, unsigned int *Triangle_NumSDsTouching, ChSystemGpu_impl::GranParamsPtr gran_params, ChSystemGpuMesh_impl::MeshParamsPtr mesh_params) |
__global__ void | storeSDsTouchedByEachTriangle (const ChSystemGpuMesh_impl::TriangleSoupPtr d_triangleSoup, const unsigned int *Triangle_NumSDsTouching, const unsigned int *TriangleSDCompositeOffsets, unsigned int *Triangle_SDsComposite, unsigned int *Triangle_TriIDsComposite, ChSystemGpu_impl::GranParamsPtr gran_params, ChSystemGpuMesh_impl::MeshParamsPtr mesh_params) |
kernel is called to populate two arrays More... | |
__device__ bool | snap_to_face (const double3 &A, const double3 &B, const double3 &C, const double3 &P, double3 &res) |
This utility function takes the location 'P' and snaps it to the closest point on the triangular face with given vertices (A, B, and C). More... | |
__device__ bool | face_sphere_cd (const double3 &A, const double3 &B, const double3 &C, const double3 &sphere_pos, const int radius, float3 &normal, float &depth, double3 &pt1) |
/brief TRIANGLE FACE - SPHERE NARROW-PHASE COLLISION DETECTION More... | |
Variables | |
const int | SAFETY_PARAM = 1000 |
Function Documentation
◆ determineCountOfSDsTouchedByEachTriangle()
__global__ void determineCountOfSDsTouchedByEachTriangle | ( | const ChSystemGpuMesh_impl::TriangleSoupPtr | d_triangleSoup, |
unsigned int * | Triangle_NumSDsTouching, | ||
ChSystemGpu_impl::GranParamsPtr | gran_params, | ||
ChSystemGpuMesh_impl::MeshParamsPtr | mesh_params | ||
) |
- Parameters
-
Triangle_NumSDsTouching number of SDs touching this Triangle
◆ face_sphere_cd()
__device__ bool face_sphere_cd | ( | const double3 & | A, |
const double3 & | B, | ||
const double3 & | C, | ||
const double3 & | sphere_pos, | ||
const int | radius, | ||
float3 & | normal, | ||
float & | depth, | ||
double3 & | pt1 | ||
) |
/brief TRIANGLE FACE - SPHERE NARROW-PHASE COLLISION DETECTION
The triangular face is defined by points A, B, C. The sequence is important as it defines the positive face via a right-hand rule. The sphere is centered at sphere_pos and has radius. The index "1" is associated with the triangle. The index "2" is associated with the sphere. The coordinates of the face and sphere are assumed to be provided in the same reference frame.
Output:
- pt1: contact point on triangle
- depth: penetration distance (a negative value means that overlap exists)
- normal: contact normal, from pt2 to pt1 A return value of "true" signals collision.
- Parameters
-
A First vertex of the triangle B Second vertex of the triangle C Third vertex of the triangle sphere_pos Location of the center of the sphere radius Sphere radius normal contact normal depth penetration pt1 contact point on triangle
◆ 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
◆ getNumberOfSpheresTouchingEachSD()
__global__ void getNumberOfSpheresTouchingEachSD | ( | ChSystemGpu_impl::GranSphereDataPtr | sphere_data, |
unsigned int | nSpheres, | ||
ChSystemGpu_impl::GranParamsPtr | gran_params | ||
) |
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 sphere touche each SD. 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_CHGPU_ID: the equivalent of a non-sphere SD ID, or a non-sphere ID
Set aside shared memory
◆ snap_to_face()
__device__ bool snap_to_face | ( | const double3 & | A, |
const double3 & | B, | ||
const double3 & | C, | ||
const double3 & | P, | ||
double3 & | res | ||
) |
This utility function takes the location 'P' and snaps it to the closest point on the triangular face with given vertices (A, B, and C).
The result is returned in 'res'. Both 'P' and 'res' are assumed to be specified in the same frame as the face vertices. This function returns 'true' if the result is on an edge of this face and 'false' if the result is inside the triangle. Code from Ericson, "real-time collision detection", 2005, pp. 141
◆ storeSDsTouchedByEachTriangle()
__global__ void storeSDsTouchedByEachTriangle | ( | const ChSystemGpuMesh_impl::TriangleSoupPtr | d_triangleSoup, |
const unsigned int * | Triangle_NumSDsTouching, | ||
const unsigned int * | TriangleSDCompositeOffsets, | ||
unsigned int * | Triangle_SDsComposite, | ||
unsigned int * | Triangle_TriIDsComposite, | ||
ChSystemGpu_impl::GranParamsPtr | gran_params, | ||
ChSystemGpuMesh_impl::MeshParamsPtr | mesh_params | ||
) |
kernel is called to populate two arrays
- Parameters
-
d_triangleSoup - the collection of triangles in this mesh soup Triangle_NumSDsTouching - number of SDs touched by each triangle TriangleSDCompositeOffsets - offsets in the array that say each SD what triangles touch it Triangle_SDsComposite - the list of SDs touched by each triangle; triangle by triangle Triangle_TriIDsComposite - array that goes hand in hand with Triangle_SDsComposite; it repeats the triangle ID for a subsequent sort by key op that is performed elsewhere gran_params mesh_params
- Returns
◆ triangle_countTouchedSDs()
|
inline |
Takes in a triangle's position in UU and finds out how many SDs it touches.
Triangle broadphase is done in float by applying the frame transform and then converting the GRF position to SU
◆ triangle_figureOutTouchedSDs()
|
inline |
Takes in a triangle's position in UU and finds out what SDs it touches.
Triangle broadphase is done in float by applying the frame transform and then converting the GRF position to SU