Functions | |
__device__ void | gpuAddHarmonicBondForce (FLOAT_ARRAY_TYPE *dForce, FLOAT_ARRAY_TYPE *dA, FLOAT_ARRAY_TYPE *dB, float *dK, float *boxsize) |
Add harmonic force to one particle. | |
__global__ void | gpuHarmonicBondForces (FLOAT_ARRAY_TYPE *dPos, FLOAT_ARRAY_TYPE *dForce, int *cnn, size_t pitch, int N, int numNN, float k, float boxsize) |
Harmonic force calculation kernel. | |
__global__ void | gpuHarmonicBondForcesSmem (FLOAT_ARRAY_TYPE *dPos, FLOAT_ARRAY_TYPE *dForce, int *cnn, size_t pitch, int N, int numNN, float k, float boxsize) |
Harmonic force calculation kernel using shared memory. |
This module contains functions to compute harmonic forces.
__device__ void gpuAddHarmonicBondForce | ( | FLOAT_ARRAY_TYPE * | dForce, | |
FLOAT_ARRAY_TYPE * | dA, | |||
FLOAT_ARRAY_TYPE * | dB, | |||
float * | dK, | |||
float * | boxsize | |||
) |
This function computes the harmonic interaction force of two particles.
[out] | dForce | Force contribution to particle a due to harmonic interaction with particle b |
[in] | dA | Position in device memory of particle a |
[in] | dB | Position in device memory of particleb |
[in] | dK | Spring constant |
[in] | boxsize | Edge length of the cubic simulation box. |
Definition at line 58 of file gpu_forces.cu.
References FLOAT_ARRAY_TYPE, and gpuDistPBC().
Referenced by gpuHarmonicBondForces(), and gpuHarmonicBondForcesSmem().
__global__ void gpuHarmonicBondForces | ( | FLOAT_ARRAY_TYPE * | dPos, | |
FLOAT_ARRAY_TYPE * | dForce, | |||
int * | cnn, | |||
size_t | pitch, | |||
int | N, | |||
int | numNN, | |||
float | k, | |||
float | boxsize | |||
) |
The harmonic force on each particle is calculated by storing the coordinates of the reference particle in a register and fetching the coordinates of the interaction partners directly from global memory.
[in] | dPos | Array of all particle positions |
[in,out] | dForce | Array of all forces |
[in] | cnn | Array holding in each column the harmonic interaction partners for one particle |
[in] | pitch | Pitch value for accessing the cnn array via CUDA in the right way. |
[in] | N | Number of particle |
[in] | numNN | Maximum number of harmonic interaction partners for each particle |
[in] | k | Spring constant |
[in] | boxsize | Edge length of the cubic simulation box |
Definition at line 422 of file gpu_forces.cu.
References FLOAT_ARRAY_TYPE, and gpuAddHarmonicBondForce().
__global__ void gpuHarmonicBondForcesSmem | ( | FLOAT_ARRAY_TYPE * | dPos, | |
FLOAT_ARRAY_TYPE * | dForce, | |||
int * | cnn, | |||
size_t | pitch, | |||
int | N, | |||
int | numNN, | |||
float | k, | |||
float | boxsize | |||
) |
The harmonic force on each particle is calculated by storing the coordinates of the reference particle in shared memory and fetching the coordinates of the interaction partners directly from global memory.
[in] | dPos | Array of all particle positions |
[in,out] | dForce | Array of all forces |
[in] | cnn | Array holding in each column the harmonic interaction partners for one particle |
[in] | pitch | Pitch value for accessing the cnn array via CUDA in the right way. |
[in] | N | Number of particle |
[in] | numNN | Maximum number of harmonic interaction partners for each particle |
[in] | k | Spring constant |
[in] | boxsize | Edge length of the cubic simulation box |
Definition at line 443 of file gpu_forces.cu.
References FLOAT_ARRAY_TYPE, and gpuAddHarmonicBondForce().