OpenMPCD
|
Contains CUDA Device code. More...
Namespaces | |
NormalMode | |
Namespace for CUDA Device functionality related to normal modes. | |
Functions | |
template<typename T > | |
__global__ void | arithmeticMean_kernel (const T *const values, const std::size_t numberOfValues, T *const output) |
Kernel to compute the arithmetic mean of the given values. More... | |
template<typename T > | |
void | arithmeticMean (const T *const values, const std::size_t numberOfValues, T *const output) |
Computes, on the CUDA Device, the arithmetic mean of the given values. More... | |
void | setLeesEdwardsSymbols (const FP shearRate, const unsigned int simBoxY) |
Sets CUDA symbols for Lees-Edwards boundary conditions. More... | |
const __device__ Vector3D< MPCParticlePositionType > | getImageUnderLeesEdwardsBoundaryConditions (const FP mpcTime, const Vector3D< MPCParticlePositionType > &position, MPCParticleVelocityType &velocityCorrection) |
Returns the image of the given particle position under Lees-Edwards boundary conditions. More... | |
__device__ bool | isInPrimarySimulationVolume (const Vector3D< MPCParticlePositionType > &position) |
Returns whether the given position lies within the primary simulation volume. More... | |
__device__ unsigned int | getCollisionCellIndex (const Vector3D< MPCParticlePositionType > &position) |
Returns the collision cell index for the given position. More... | |
__global__ void | sortParticlesIntoCollisionCellsLeesEdwards (const unsigned int workUnitOffset, const unsigned int particleCount, const MPCParticlePositionType *const gridShift_, const FP mpcTime, const MPCParticlePositionType *const positions, MPCParticleVelocityType *const velocities, MPCParticleVelocityType *const velocityCorrections, unsigned int *const collisionCellIndices) |
Sorts the MPC particles into the collision cells, temporarily applying Lees-Edwards boundary conditions. More... | |
__device__ void | sortIntoCollisionCellsLeesEdwards (const unsigned int particleID, const MPCParticlePositionType *const gridShift_, const FP mpcTime, const MPCParticlePositionType *const positions, MPCParticleVelocityType *const velocities, MPCParticleVelocityType *const velocityCorrections, unsigned int *const collisionCellIndices) |
Sorts the given particle into the collision cells, temporarily applying Lees-Edwards boundary conditions. More... | |
__global__ void | collisionCellContributions (const unsigned int workUnitOffset, const unsigned int particleCount, const MPCParticleVelocityType *const velocities, const unsigned int *const collisionCellIndices, MPCParticleVelocityType *const collisionCellMomenta, FP *const collisionCellMasses, const FP particleMass) |
Computes the collision cell mass and momentum contributions by the given particles. More... | |
__device__ Vector3D< MPCParticleVelocityType > | getCollisionCellCenterOfMassVelocity (const unsigned int collisionCellIndex, const MPCParticleVelocityType *const collisionCellMomenta, const FP *const collisionCellMasses) |
Returns the center-of-mass velocity of a collision cell. More... | |
__global__ void | resetCollisionCellData (const unsigned int workUnitOffset, const unsigned int collisionCellCount, unsigned int *const collisionCellParticleCounts, MPCParticleVelocityType *const collisionCellMomenta, FP *const collisionCellFrameInternalKineticEnergies, FP *const collisionCellMasses) |
Resets the collision cell data buffers. More... | |
__global__ void | generateGridShiftVector (MPCParticlePositionType *const gridShift, const FP gridShiftScale, GPURNG *const rngs) |
Generates a new random grid shift vector. More... | |
__global__ void | generateCollisionCellRotationAxes (const unsigned int workUnitOffset, const unsigned int collisionCellCount, FP *const collisionCellRotationAxes, GPURNG *const rngs) |
Generates new rotation axes for the collision cells. More... | |
__global__ void | collisionCellStochasticRotationStep1 (const unsigned int workUnitOffset, const unsigned int particleCount, MPCParticleVelocityType *const velocities, const FP particleMass, const unsigned int *const collisionCellIndices, const MPCParticleVelocityType *const collisionCellMomenta, const FP *const collisionCellMasses, const MPCParticlePositionType *const collisionCellRotationAxes, FP *const collisionCellFrameInternalKineticEnergies, unsigned int *const collisionCellParticleCounts) |
Applies the first step of the SRD rotation to the given particles. More... | |
__global__ void | generateCollisionCellMBSFactors (const unsigned int workUnitOffset, const unsigned int collisionCellCount, const FP *const collisionCellFrameInternalKineticEnergies, unsigned int *const collisionCellParticleCounts, FP *const collisionCellRelativeVelocityScalings, const FP bulkThermostatTargetkT, GPURNG *const rngs) |
Generates Maxwell-Boltzmann-Scaling factors for the collision cells. More... | |
__global__ void | collisionCellStochasticRotationStep2 (const unsigned int workUnitOffset, const unsigned int particleCount, MPCParticleVelocityType *const velocities, const unsigned int *const collisionCellIndices, const MPCParticleVelocityType *const collisionCellMomenta, const FP *const collisionCellMasses, const FP *const collisionCellVelocityScalings) |
Applies the second step of the SRD rotation to the given particles. More... | |
__global__ void | undoLeesEdwardsVelocityCorrections (const unsigned int workUnitOffset, const unsigned int particleCount, MPCParticleVelocityType *const velocities, const MPCParticleVelocityType *const velocityCorrections) |
Undoes the velocity corrections applied by sortIntoCollisionCellsLeesEdwards . More... | |
__global__ void | constructGPURNGs (const std::size_t count, GPURNG *const location, const unsigned long long seed) |
Sets up instances of GPURNG in the specified memory location. More... | |
__global__ void | destroyGPURNGs (const std::size_t count, GPURNG *const location) |
Destroys instances of GPURNG in the specified memory location. More... | |
void | setSimulationBoxSizeSymbols (const unsigned int simBoxSizeX, const unsigned int simBoxSizeY, const unsigned int simBoxSizeZ) |
Sets the symbols mpcSimulationBoxSizeX , mpcSimulationBoxSizeY , mpcSimulationBoxSizeZ , and collisionCellCount in OpenMPCD::CUDA::DeviceCode . More... | |
void | setMPCStreamingTimestep (const FP timestep) |
Sets the symbol OpenMPCD::CUDA::DeviceCode::streamingTimestep . More... | |
void | setSRDCollisionAngleSymbol (const FP angle) |
Sets the symbol OpenMPCD::CUDA::DeviceCode::srdCollisionAngle . More... | |
OPENMPCD_CUDA_DEVICE void | atomicAdd (double *const target, const double increment) |
Atomically adds increment to target. More... | |
OPENMPCD_CUDA_DEVICE void | atomicAdd (float *const target, const float increment) |
Atomically adds increment to target. More... | |
template<typename B > | |
OPENMPCD_CUDA_HOST_AND_DEVICE boost::enable_if< boost::is_integral< B >, double >::type | pow (const B base, const double exponent) |
The power function. More... | |
OPENMPCD_CUDA_HOST_AND_DEVICE double | pow (const double base, const double exponent) |
The power function. More... | |
__device__ MPCParticlePositionType | velocityVerletStep1 (const MPCParticlePositionType position, const MPCParticleVelocityType velocity, const FP acceleration, const FP timestep) |
Performs the first step in the velocity-Verlet algorithm. More... | |
__device__ MPCParticleVelocityType | velocityVerletStep2 (const MPCParticleVelocityType velocity, const FP oldAcceleration, const FP newAcceleration, const FP timestep) |
Performs the second step in the velocity-Verlet algorithm. More... | |
__device__ void | velocityVerletStep1 (RemotelyStoredVector< MPCParticlePositionType > *const position, const RemotelyStoredVector< MPCParticleVelocityType > velocity, const Vector3D< FP > acceleration, const FP timestep) |
Performs the first step in the velocity-Verlet algorithm. More... | |
__device__ void | velocityVerletStep2 (RemotelyStoredVector< MPCParticleVelocityType > *const velocity, const Vector3D< FP > oldAcceleration, const Vector3D< FP > newAcceleration, const FP timestep) |
Performs the second step in the velocity-Verlet algorithm. More... | |
Variables | |
__constant__ FP | leesEdwardsRelativeLayerVelocity |
The relative x -velocity of y -adjacent layers in Lees-Edwards boundary conditions. More... | |
__constant__ unsigned int | mpcSimulationBoxSizeX |
The size of the primary simulation box along the x direction. More... | |
__constant__ unsigned int | mpcSimulationBoxSizeY |
The size of the primary simulation box along the y direction. More... | |
__constant__ unsigned int | mpcSimulationBoxSizeZ |
The size of the primary simulation box along the z direction. More... | |
__constant__ unsigned int | collisionCellCount |
The number of collision cells in the system. More... | |
__constant__ FP | streamingTimestep |
The MPC streaming time-step. More... | |
__constant__ FP | srdCollisionAngle |
The collision angle for SRD collisions. More... | |
void OpenMPCD::CUDA::DeviceCode::arithmeticMean | ( | const T *const | values, |
const std::size_t | numberOfValues, | ||
T *const | output | ||
) |
Computes, on the CUDA Device, the arithmetic mean of the given values.
Each value will enter the average with equal weight.
OpenMPCD::NULLPointerException | Throws if values is nullptr . |
OpenMPCD::InvalidArgumentException | Throws if numberOfValues == 0 . |
OpenMPCD::NULLPointerException | Throws if output is nullptr . |
OpenMPCD::InvalidArgumentException | If OPENMPCD_DEBUG is defined, throws if values is not a Device pointer. |
OpenMPCD::InvalidArgumentException | If OPENMPCD_DEBUG is defined, throws if output is not a Device pointer. |
T | The type of values. |
[in] | values | The values to take the average of, stored on the CUDA Device. In total, values must hold at least numberOfValues * sizeof(T) bytes. |
[in] | numberOfValues | The number of values to treat. |
[out] | output | The CUDA Device pointer to save the result to. |
Definition at line 33 of file ImplementationDetails/Average.hpp.
__global__ void OpenMPCD::CUDA::DeviceCode::arithmeticMean_kernel | ( | const T *const | values, |
const std::size_t | numberOfValues, | ||
T *const | output | ||
) |
Kernel to compute the arithmetic mean of the given values.
Each value will enter the average with equal weight.
T | The type of values. |
[in] | values | The values to take the average of, stored on the CUDA Device. In total, values must hold at least numberOfValues * sizeof(T) bytes. |
[in] | numberOfValues | The number of values to treat. |
[out] | output | The CUDA Device pointer to save the result to. |
Definition at line 22 of file ImplementationDetails/Average.hpp.
OPENMPCD_CUDA_DEVICE void OpenMPCD::CUDA::DeviceCode::atomicAdd | ( | double *const | target, |
const double | increment | ||
) |
Atomically adds increment to target.
[in] | target | The address of the value to increment. |
[in] | increment | The value to add to target. |
OPENMPCD_CUDA_DEVICE void OpenMPCD::CUDA::DeviceCode::atomicAdd | ( | float *const | target, |
const float | increment | ||
) |
Atomically adds increment to target.
[in] | target | The address of the value to increment. |
[in] | increment | The value to add to target. |
__global__ void OpenMPCD::CUDA::DeviceCode::collisionCellContributions | ( | const unsigned int | workUnitOffset, |
const unsigned int | particleCount, | ||
const MPCParticleVelocityType *const | velocities, | ||
const unsigned int *const | collisionCellIndices, | ||
MPCParticleVelocityType *const | collisionCellMomenta, | ||
FP *const | collisionCellMasses, | ||
const FP | particleMass | ||
) |
Computes the collision cell mass and momentum contributions by the given particles.
[in] | workUnitOffset | The number of particles to skip. |
[in] | particleCount | The number of particles there are. |
[in] | velocities | The particle velocities. |
[in] | collisionCellIndices | The collision cell indices for the particles. |
[in,out] | collisionCellMomenta | The collision cell momenta. This array is expected to be long enough to store 3 coordinates for each collision cell, and is assumed to have each element set to 0 prior to calling this function. |
[in,out] | collisionCellMasses | The masses of the collision cells. This array is expected to be long enough to store 1 element for each collision cell, and is assumed to have each element set to 0 prior to calling this function. |
[in] | particleMass | The mass of any one particle. |
__global__ void OpenMPCD::CUDA::DeviceCode::collisionCellStochasticRotationStep1 | ( | const unsigned int | workUnitOffset, |
const unsigned int | particleCount, | ||
MPCParticleVelocityType *const | velocities, | ||
const FP | particleMass, | ||
const unsigned int *const | collisionCellIndices, | ||
const MPCParticleVelocityType *const | collisionCellMomenta, | ||
const FP *const | collisionCellMasses, | ||
const MPCParticlePositionType *const | collisionCellRotationAxes, | ||
FP *const | collisionCellFrameInternalKineticEnergies, | ||
unsigned int *const | collisionCellParticleCounts | ||
) |
Applies the first step of the SRD rotation to the given particles.
This function requires that OpenMPCD::CUDA::DeviceCode::setSRDCollisionAngleSymbol
has been called before.
[in] | workUnitOffset | The number of particles to skip. |
[in] | particleCount | The number of particles there are. |
[in,out] | velocities | The particle velocities as input; after the function call, the values in this array will contain the velocities in the collision cell's center-of-mass frame. |
[in] | particleMass | The mass of any one particle. |
[in] | collisionCellIndices | An array which holds, for each particle, the index of the collision cell it is currently in. |
[in] | collisionCellMomenta | An array holding, for each collision cell, its three cartesian momentum coordinates. |
[in] | collisionCellMasses | An array holding, for each collision cell, its total mass content. |
[in] | collisionCellRotationAxes | An array holding, for each collision cell, the random rotation axis it is to use during the collision. |
[out] | collisionCellFrameInternalKineticEnergies | An array, long enough to hold an element for each collision cell. The respective elements are set to the sum of the kinetic energies of the particles in that collision cell, as measured in that collision cell's center-of-mass frame. The array is not set to 0 by this function prior to execution of the algorithm. |
[out] | collisionCellParticleCounts | An array, long enough to hold an element for each collision cell. The respective elements are set to the number of particles in that collision cell. The array is assumed to contain only values 0 prior to calling this function. |
__global__ void OpenMPCD::CUDA::DeviceCode::collisionCellStochasticRotationStep2 | ( | const unsigned int | workUnitOffset, |
const unsigned int | particleCount, | ||
MPCParticleVelocityType *const | velocities, | ||
const unsigned int *const | collisionCellIndices, | ||
const MPCParticleVelocityType *const | collisionCellMomenta, | ||
const FP *const | collisionCellMasses, | ||
const FP *const | collisionCellVelocityScalings | ||
) |
Applies the second step of the SRD rotation to the given particles.
This kernel scales each particle's velocity with the corresponding collision cell's scaling factor (stored in collisionCellVelocityScalings
), and then adds the corresponding collision cell's center-of-mass velocity.
[in] | workUnitOffset | The number of particles to skip. |
[in] | particleCount | The number of particles there are. |
[in,out] | velocities | The new particle velocities as output; as input, the values in this array are to contain the rotated velocities in the collision cell's center-of-mass frame. |
[in] | collisionCellIndices | An array which holds, for each particle, the index of the collision cell it is currently in. |
[in] | collisionCellMomenta | An array holding, for each collision cell, its three Cartesian momentum coordinates. |
[in] | collisionCellMasses | An array holding, for each collision cell, its total mass content. |
[in] | collisionCellVelocityScalings | For each collision cell, stores the factor by which to scale the relative velocities. |
__global__ void OpenMPCD::CUDA::DeviceCode::constructGPURNGs | ( | const std::size_t | count, |
GPURNG *const | location, | ||
const unsigned long long | seed | ||
) |
Sets up instances of GPURNG
in the specified memory location.
The individual RNGs will have subsequence numbers (or, alternatively, seeds) set such that the individual instances generate independent streams of random numbers.
Only the first thread in this kernel does any work.
[in] | count | The number of instances to construct. |
[out] | location | The location in memory to store the instances in; memory must have been allocated at this location prior to calling this function. |
[in] | seed | The seed to use for the RNGs. |
__global__ void OpenMPCD::CUDA::DeviceCode::destroyGPURNGs | ( | const std::size_t | count, |
GPURNG *const | location | ||
) |
Destroys instances of GPURNG
in the specified memory location.
Only the first thread in this kernel does any work.
[in] | count | The number of instances to construct. |
[out] | location | The location in memory where there are instances; memory must be freed manually after this function call. |
__global__ void OpenMPCD::CUDA::DeviceCode::generateCollisionCellMBSFactors | ( | const unsigned int | workUnitOffset, |
const unsigned int | collisionCellCount, | ||
const FP *const | collisionCellFrameInternalKineticEnergies, | ||
unsigned int *const | collisionCellParticleCounts, | ||
FP *const | collisionCellRelativeVelocityScalings, | ||
const FP | bulkThermostatTargetkT, | ||
GPURNG *const | rngs | ||
) |
Generates Maxwell-Boltzmann-Scaling factors for the collision cells.
[in] | workUnitOffset | The number of particles to skip. |
[in] | collisionCellCount | The number of collision cells there are. |
[in] | collisionCellFrameInternalKineticEnergies | Stores, for each collision cell, the sum of the kinetic energies of the particles in that collision cell, as measured in that collision cell's center-of-mass frame |
[in] | collisionCellParticleCounts | The number of particles in the collision cells. |
[out] | collisionCellRelativeVelocityScalings | For each collision cell, stores the factor by which to scale the relative velocities. |
[in] | bulkThermostatTargetkT | The target temperature for the thermostat. |
[in,out] | rngs | The random number generators, of which there must be at least collisionCellCount . |
__global__ void OpenMPCD::CUDA::DeviceCode::generateCollisionCellRotationAxes | ( | const unsigned int | workUnitOffset, |
const unsigned int | collisionCellCount, | ||
FP *const | collisionCellRotationAxes, | ||
GPURNG *const | rngs | ||
) |
Generates new rotation axes for the collision cells.
[in] | workUnitOffset | The number of particles to skip. |
[in] | collisionCellCount | The number of collision cells there are. |
[out] | collisionCellRotationAxes | For each collision cell, stores the Cartesian coordinates of the unit-length axis of rotation. |
[in,out] | rngs | The random number generators, of which there must be at least collisionCellCount . |
__global__ void OpenMPCD::CUDA::DeviceCode::generateGridShiftVector | ( | MPCParticlePositionType *const | gridShift, |
const FP | gridShiftScale, | ||
GPURNG *const | rngs | ||
) |
Generates a new random grid shift vector.
This function will draw three Cartesian coordinates from the uniform distribution over \( \left( 0, 1 \right] \), multiply the results by gridShiftScale
, and store them in the three elements pointed at by gridShift
.
This kernel must be called with one block and one thread only.
[out] | gridShift | Where to store the three grid shift coordinates to. |
[in] | gridShiftScale | The scale for the grid shift vector coordinates. |
[in,out] | rngs | Pointer to at least one random number generator. |
__device__ Vector3D<MPCParticleVelocityType> OpenMPCD::CUDA::DeviceCode::getCollisionCellCenterOfMassVelocity | ( | const unsigned int | collisionCellIndex, |
const MPCParticleVelocityType *const | collisionCellMomenta, | ||
const FP *const | collisionCellMasses | ||
) |
Returns the center-of-mass velocity of a collision cell.
[in] | collisionCellIndex | The index of the collision cell. |
[in] | collisionCellMomenta | The momenta of the collision cells. |
[in] | collisionCellMasses | The masses of the collision cells. |
Definition at line 113 of file Simulation.cpp.
__device__ unsigned int OpenMPCD::CUDA::DeviceCode::getCollisionCellIndex | ( | const Vector3D< MPCParticlePositionType > & | position | ) |
Returns the collision cell index for the given position.
The collision cell index is determined as follows: The three cartesian coordinates of the given particle position are rounded down towards the nearest integer. This triple of numbers then determine the coordinates of the collision cell the particle lies in. If they are called cellX
, cellY
, and cellZ
, respectively, what is returned is cellX + cellY * boxSizeX
+ cellZ * boxSizeX * boxSizeY, where
boxSizeX and
boxSizeYare the number of collision cells in the primary simulation volume along the
xand
y` directions, respectively.
This function requires that OpenMPCD::CUDA::DeviceCode::setSimulationBoxSizeSymbols
has been called before.
[in] | position | The position to consider, which must lie in the primary simulation volume, i.e. the x , y , and z coordinates must be non-negative and smaller than the simulation box size along the respective direction. |
const __device__ Vector3D<MPCParticlePositionType> OpenMPCD::CUDA::DeviceCode::getImageUnderLeesEdwardsBoundaryConditions | ( | const FP | mpcTime, |
const Vector3D< MPCParticlePositionType > & | position, | ||
MPCParticleVelocityType & | velocityCorrection | ||
) |
Returns the image of the given particle position under Lees-Edwards boundary conditions.
This function assumes that setLeesEdwardsSymbols
and setSimulationBoxSizeSymbols
have been called before.
[in] | mpcTime | The simulation time for the MPC steps. |
[in] | position | The particle position. |
[out] | velocityCorrection | Set to the velocity along the x direction that needs to be added to the particle's velocity. |
__device__ bool OpenMPCD::CUDA::DeviceCode::isInPrimarySimulationVolume | ( | const Vector3D< MPCParticlePositionType > & | position | ) |
Returns whether the given position lies within the primary simulation volume.
The primary simulation volume is defined as
\[ \left[0, L_x\right) \times \left[0, L_y\right) \times \left[0, L_z\right) \]
where \( L_x \) is the value of mpcSimulationBoxSizeX
, i.e. the size of the primary simulation volume along the \( x \) axis, and analogously for \( L_y \) and mpcSimulationBoxSizeY
, and \( L_z \) and mpcSimulationBoxSizeZ
.
This function requires that OpenMPCD::CUDA::DeviceCode::setSimulationBoxSizeSymbols
has been called.
[in] | position | The position to check. |
OPENMPCD_CUDA_HOST_AND_DEVICE boost::enable_if<boost::is_integral<B>, double>::type OpenMPCD::CUDA::DeviceCode::pow | ( | const B | base, |
const double | exponent | ||
) |
The power function.
B | The type of the base argument, which must be an integral type. |
[in] | base | The base. |
[in] | exponent | The exponent. |
Definition at line 50 of file Utilities.hpp.
|
inline |
The power function.
[in] | base | The base. |
[in] | exponent | The exponent. |
Definition at line 62 of file Utilities.hpp.
__global__ void OpenMPCD::CUDA::DeviceCode::resetCollisionCellData | ( | const unsigned int | workUnitOffset, |
const unsigned int | collisionCellCount, | ||
unsigned int *const | collisionCellParticleCounts, | ||
MPCParticleVelocityType *const | collisionCellMomenta, | ||
FP *const | collisionCellFrameInternalKineticEnergies, | ||
FP *const | collisionCellMasses | ||
) |
Resets the collision cell data buffers.
This kernel sets the data pointed to by the arguments provided to 0
. One needs to spawn at least collisionCellCount
threads (possibly using workUnitOffset
, if one grid does not fit all threads) for this operation to complete semantically.
Each pointer argument is required to be non-nullptr
, and point to at least collisionCellCount
elements; collisionCellMomenta
needs to point to at least 3 * collisionCellCount
elements instead.
[in] | workUnitOffset | The number of collision cells to skip. |
[in] | collisionCellCount | The number of collision cells there are. |
[out] | collisionCellParticleCounts | Pointer to the memory where the collision cell particle counts are stored. |
[out] | collisionCellMomenta | Pointer to the memory where the collision cell momenta are stored (three per collision cell). |
[out] | collisionCellFrameInternalKineticEnergies | Pointer to the memory where the collision cell kinetic energies (in the internal frame) are stored. |
[out] | collisionCellMasses | Pointer to the memory where the collision cell masses are stored. |
void OpenMPCD::CUDA::DeviceCode::setLeesEdwardsSymbols | ( | const FP | shearRate, |
const unsigned int | simBoxY | ||
) |
Sets CUDA symbols for Lees-Edwards boundary conditions.
This function stores \( \Delta v_x = \dot{\gamma} L_y \) in leesEdwardsRelativeLayerVelocity
.
[in] | shearRate | The shear rate \( \dot{\gamma} \). |
[in] | simBoxY | The size \( L_y \) of the primary simulation box along the y direction. |
void OpenMPCD::CUDA::DeviceCode::setMPCStreamingTimestep | ( | const FP | timestep | ) |
Sets the symbol OpenMPCD::CUDA::DeviceCode::streamingTimestep
.
OpenMPCD::InvalidArgumentException | If OPENMPCD_DEBUG is defined, throws if timestep == 0 . |
[in] | timestep | The streaming timestep for each MPC fluid streaming step. |
void OpenMPCD::CUDA::DeviceCode::setSimulationBoxSizeSymbols | ( | const unsigned int | simBoxSizeX, |
const unsigned int | simBoxSizeY, | ||
const unsigned int | simBoxSizeZ | ||
) |
Sets the symbols mpcSimulationBoxSizeX
, mpcSimulationBoxSizeY
, mpcSimulationBoxSizeZ
, and collisionCellCount
in OpenMPCD::CUDA::DeviceCode
.
collisionCellCount
will be set to the product of the three arguments supplied.
OpenMPCD::InvalidArgumentException | If OPENMPCD_DEBUG is defined, throws if either argument is 0 . |
[in] | simBoxSizeX | The value for mpcSimulationBoxSizeX . |
[in] | simBoxSizeY | The value for mpcSimulationBoxSizeY . |
[in] | simBoxSizeZ | The value for mpcSimulationBoxSizeZ . |
void OpenMPCD::CUDA::DeviceCode::setSRDCollisionAngleSymbol | ( | const FP | angle | ) |
Sets the symbol OpenMPCD::CUDA::DeviceCode::srdCollisionAngle
.
[in] | angle | The SRD collision angle. |
__device__ void OpenMPCD::CUDA::DeviceCode::sortIntoCollisionCellsLeesEdwards | ( | const unsigned int | particleID, |
const MPCParticlePositionType *const | gridShift_, | ||
const FP | mpcTime, | ||
const MPCParticlePositionType *const | positions, | ||
MPCParticleVelocityType *const | velocities, | ||
MPCParticleVelocityType *const | velocityCorrections, | ||
unsigned int *const | collisionCellIndices | ||
) |
Sorts the given particle into the collision cells, temporarily applying Lees-Edwards boundary conditions.
This function requires that OpenMPCD::CUDA::DeviceCode::setSimulationBoxSizeSymbols
and OpenMPCD::CUDA::DeviceCode::setLeesEdwardsSymbols
have been called before.
[in] | particleID | The particle ID. |
[in] | gridShift_ | The grid shift vector to temporarily add to the position. |
[in] | mpcTime | The MPC time that has passed since the start of the simulation. |
[in] | positions | The particle positions. |
[in,out] | velocities | The particle velocities, which may be changed due to Lees-Edwards boundary conditions. |
[out] | velocityCorrections | The velocity corrections applied due to Lees-Edwards boundary conditions. |
[out] | collisionCellIndices | The collision cell indices for the particles. |
__global__ void OpenMPCD::CUDA::DeviceCode::sortParticlesIntoCollisionCellsLeesEdwards | ( | const unsigned int | workUnitOffset, |
const unsigned int | particleCount, | ||
const MPCParticlePositionType *const | gridShift_, | ||
const FP | mpcTime, | ||
const MPCParticlePositionType *const | positions, | ||
MPCParticleVelocityType *const | velocities, | ||
MPCParticleVelocityType *const | velocityCorrections, | ||
unsigned int *const | collisionCellIndices | ||
) |
Sorts the MPC particles into the collision cells, temporarily applying Lees-Edwards boundary conditions.
This function requires that OpenMPCD::CUDA::DeviceCode::setSimulationBoxSizeSymbols
and OpenMPCD::CUDA::DeviceCode::setLeesEdwardsSymbols
have been called before.
[in] | workUnitOffset | The number of particles to skip. |
[in] | particleCount | The number of particles there are. |
[in] | gridShift_ | The grid shift vector to temporarily add to all positions. |
[in] | mpcTime | The MPC time that has passed since the start of the simulation. |
[in] | positions | The MPC fluid particle positions. |
[in,out] | velocities | The MPC fluid particle velocities, which may be changed due to Lees-Edwards boundary conditions. |
[out] | velocityCorrections | The velocity corrections applied due to Lees-Edwards boundary conditions, which can be undone by undoLeesEdwardsVelocityCorrections ; the buffer has to be at least particleCount elements long. |
[out] | collisionCellIndices | The collision cell indices for the particles. |
__global__ void OpenMPCD::CUDA::DeviceCode::undoLeesEdwardsVelocityCorrections | ( | const unsigned int | workUnitOffset, |
const unsigned int | particleCount, | ||
MPCParticleVelocityType *const | velocities, | ||
const MPCParticleVelocityType *const | velocityCorrections | ||
) |
Undoes the velocity corrections applied by sortIntoCollisionCellsLeesEdwards
.
[in] | workUnitOffset | The number of particles to skip. |
[in] | particleCount | The number of particles there are. |
[in,out] | velocities | The particle velocities. |
[in] | velocityCorrections | The velocity corrections as returned by sortIntoCollisionCellsLeesEdwards . |
__device__ MPCParticlePositionType OpenMPCD::CUDA::DeviceCode::velocityVerletStep1 | ( | const MPCParticlePositionType | position, |
const MPCParticleVelocityType | velocity, | ||
const FP | acceleration, | ||
const FP | timestep | ||
) |
Performs the first step in the velocity-Verlet algorithm.
This function returns the new position \(x(t+h)\) according to \(x(t+h) = x(t) + v(t) h + \frac{1}{2} a h^2\) where \(v\) is the velocity, \(a\) the acceleration and \(h\) the timestep.
[in] | position | The position. |
[in] | velocity | The velocity. |
[in] | acceleration | The acceleration. |
[in] | timestep | The timestep. |
__device__ void OpenMPCD::CUDA::DeviceCode::velocityVerletStep1 | ( | RemotelyStoredVector< MPCParticlePositionType > *const | position, |
const RemotelyStoredVector< MPCParticleVelocityType > | velocity, | ||
const Vector3D< FP > | acceleration, | ||
const FP | timestep | ||
) |
Performs the first step in the velocity-Verlet algorithm.
This function updates the position \(x\) according to \(x(t+h) = x(t) + v(t) h + \frac{1}{2} a h^2\) where \(v\) is the velocity, \(a\) the acceleration and \(h\) the timestep.
[in,out] | position | The position. |
[in] | velocity | The velocity. |
[in] | acceleration | The acceleration. |
[in] | timestep | The timestep. |
__device__ MPCParticleVelocityType OpenMPCD::CUDA::DeviceCode::velocityVerletStep2 | ( | const MPCParticleVelocityType | velocity, |
const FP | oldAcceleration, | ||
const FP | newAcceleration, | ||
const FP | timestep | ||
) |
Performs the second step in the velocity-Verlet algorithm.
This function returns the new velocity \(y(t+h)\) according to \(v(t+h) = v(t) + \frac{h}{2}\left(a(t) + a(t+h)\right)\) where \(x\) is the position, \(a\) the acceleration and \(h\) the timestep.
[in] | velocity | The velocity. |
[in] | oldAcceleration | The acceleration at the time \(t\). |
[in] | newAcceleration | The acceleration at the time \(t+h\). |
[in] | timestep | The timestep. |
__device__ void OpenMPCD::CUDA::DeviceCode::velocityVerletStep2 | ( | RemotelyStoredVector< MPCParticleVelocityType > *const | velocity, |
const Vector3D< FP > | oldAcceleration, | ||
const Vector3D< FP > | newAcceleration, | ||
const FP | timestep | ||
) |
Performs the second step in the velocity-Verlet algorithm.
This function updates the velocity \(v\) according to \(v(t+h) = v(t) + \frac{h}{2}\left(a(t) + a(t+h)\right)\) where \(x\) is the position, \(a\) the acceleration and \(h\) the timestep.
[in,out] | velocity | The velocity. |
[in] | oldAcceleration | The acceleration at the time \(t\). |
[in] | newAcceleration | The acceleration at the time \(t+h\). |
[in] | timestep | The timestep. |
__constant__ unsigned int OpenMPCD::CUDA::DeviceCode::collisionCellCount |
The number of collision cells in the system.
__constant__ FP OpenMPCD::CUDA::DeviceCode::leesEdwardsRelativeLayerVelocity |
The relative x
-velocity of y
-adjacent layers in Lees-Edwards boundary conditions.
__constant__ unsigned int OpenMPCD::CUDA::DeviceCode::mpcSimulationBoxSizeX |
The size of the primary simulation box along the x
direction.
__constant__ unsigned int OpenMPCD::CUDA::DeviceCode::mpcSimulationBoxSizeY |
The size of the primary simulation box along the y
direction.
__constant__ unsigned int OpenMPCD::CUDA::DeviceCode::mpcSimulationBoxSizeZ |
The size of the primary simulation box along the z
direction.
__constant__ FP OpenMPCD::CUDA::DeviceCode::srdCollisionAngle |
The collision angle for SRD collisions.
__constant__ FP OpenMPCD::CUDA::DeviceCode::streamingTimestep |
The MPC streaming time-step.