Functions | |
float | getHarmonicPotentialEnergy (float *dEpot, FLOAT_ARRAY_TYPE *dPos, float *sum, int *cnn, size_t &pitch, int &N, int &numNN, float &k, float &boxSize, dim3 &dimGrid, dim3 &dimBlock, int &maxThreads, int &maxBlocks) |
Calculates total harmonic potential energy. | |
float | getKineticEnergy (FLOAT_ARRAY_TYPE *dVel, int &N, float *sumSq, int &maxThreads, int &maxBlocks) |
Calculates total kinetic energy. | |
float | getLJCutPotentialEnergy (FLOAT_ARRAY_TYPE *dPos, float *dEpot, float *sum, int &N, float &boxSize, float &rCutSq, dim3 &dimGrid, dim3 &dimBlock, int &maxThreads, int &maxBlocks) |
Calculates the total Lennard-Jones energy of the system. | |
__global__ void | gpuGetHarmonicPotentialEnergyKernel (float *dEpot, FLOAT_ARRAY_TYPE *dPos, int *cnn, size_t pitch, int N, int numNN, float k, float boxSize) |
Kernel for harmonic potential. | |
__global__ void | gpuLJCutPotentialEnergyKernel (FLOAT_ARRAY_TYPE *dPos, float *dEpot, int N, float boxSize, float rCutSq) |
Calculates Lennard-Jones potential within a certrain cutoff. |
This module contains basic measurements of energies
float getHarmonicPotentialEnergy | ( | float * | dEpot, | |
FLOAT_ARRAY_TYPE * | dPos, | |||
float * | sum, | |||
int * | cnn, | |||
size_t & | pitch, | |||
int & | N, | |||
int & | numNN, | |||
float & | k, | |||
float & | boxSize, | |||
dim3 & | dimGrid, | |||
dim3 & | dimBlock, | |||
int & | maxThreads, | |||
int & | maxBlocks | |||
) |
Calculates total harmonic potential energy
[out] | dEpot | Array on the device to hold the total potential energy for each particle |
[in] | dPos | Particle position array |
[in] | sum | Summation variable on the device used for intermediate data during reduction, content will be overwritten |
[in] | cnn | Array holding the indices of harmonic linked particles |
[in] | pitch | Memory alignement pitch on the device |
[in] | N | Number of particle |
[in] | numNN | maximum number of connected neighbors |
[in] | k | Spring constant |
[in] | boxSize | Edge length of the cubic simulation box |
[in] | dimGrid | Grid dimension |
[in] | dimBlock | Block dimension |
[in] | maxThreads | Maximum number of allowed threads |
[in] | maxBlocks | Maximum number of allowed blocks |
Definition at line 269 of file gpu_tools.cu.
References checkCudaError(), and sumReduction().
float getKineticEnergy | ( | FLOAT_ARRAY_TYPE * | dVel, | |
int & | N, | |||
float * | sumSq, | |||
int & | maxThreads, | |||
int & | maxBlocks | |||
) |
Calculates total kinetic energy , where
due to Lennard-Jones units.
[in] | dVel | Velocity array |
[in] | N | Number of particles |
[in] | sumSq | Helper variable, content will be overwritten |
[in] | maxThreads | Maximum number of allowed threads |
[in] | maxBlocks | Maximum number of allowed blocks |
Definition at line 255 of file gpu_tools.cu.
References SqVecSumReduction().
float getLJCutPotentialEnergy | ( | FLOAT_ARRAY_TYPE * | dPos, | |
float * | dEpot, | |||
float * | sum, | |||
int & | N, | |||
float & | boxSize, | |||
float & | rCutSq, | |||
dim3 & | dimGrid, | |||
dim3 & | dimBlock, | |||
int & | maxThreads, | |||
int & | maxBlocks | |||
) |
Calculates the total Lennard-Jones energy of the system: , where due to Lennard-Jones units
[in] | dPos | Particle position array |
[in] | dEpot | Array on the device to hold the total potential energy for each particle |
[in] | sum | Summation variable on the device used for intermediate data during reduction, content will be overwritten |
[in] | N | Number of particle |
[in] | boxSize | Edge length of the cubic simulation box |
[in] | rCutSq | Squared Lennard-Jones cutoff radius |
[in] | dimGrid | Grid dimension |
[in] | dimBlock | Block dimension |
[in] | maxThreads | Maximum number of allowed threads |
[in] | maxBlocks | Maximum number of allowed blocks |
Definition at line 303 of file gpu_tools.cu.
References checkCudaError(), and sumReduction().
__global__ void gpuGetHarmonicPotentialEnergyKernel | ( | float * | dEpot, | |
FLOAT_ARRAY_TYPE * | dPos, | |||
int * | cnn, | |||
size_t | pitch, | |||
int | N, | |||
int | numNN, | |||
float | k, | |||
float | boxSize | |||
) |
Kernel function called by getHarmonicPotentialEnergy to calculate the harmonic potential.
[out] | dEpot | Array on the device to hold the total potential energy for each particle |
[in] | dPos | Particle position array |
[in] | cnn | Array holding the indices of harmonic linked particles |
[in] | pitch | Memory alignement pitch on the device |
[in] | N | Number of particle |
[in] | numNN | maximum number of connected neighbors |
[in] | k | Spring constant |
[in] | boxSize | Edge length of the cubic simulation box |
Definition at line 280 of file gpu_tools.cu.
References FLOAT_ARRAY_TYPE, and gpuDistPBC().
__global__ void gpuLJCutPotentialEnergyKernel | ( | FLOAT_ARRAY_TYPE * | dPos, | |
float * | dEpot, | |||
int | N, | |||
float | boxSize, | |||
float | rCutSq | |||
) |
Calculates the total Lennard-Jones energy of the system like getLJCutPotentialEnergy, but sums only contributions with .
[in] | dPos | Particle position array |
[out] | dEpot | Array on the device to hold the total potential energy for each particle |
[in] | N | Number of particle |
[in] | boxSize | Edge length of the cubic simulation box |
[in] | rCutSq | Squared Lennard-Jones cutoff radius |
Definition at line 314 of file gpu_tools.cu.
References FLOAT_ARRAY_TYPE, and gpuDistPBC().