OpenMPCD
Namespaces | Functions | Variables
OpenMPCD::CUDA::DeviceCode Namespace Reference

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< MPCParticlePositionTypegetImageUnderLeesEdwardsBoundaryConditions (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< MPCParticleVelocityTypegetCollisionCellCenterOfMassVelocity (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...
 

Detailed Description

Contains CUDA Device code.

Function Documentation

◆ arithmeticMean()

template<typename T >
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.

Exceptions
OpenMPCD::NULLPointerExceptionThrows if values is nullptr.
OpenMPCD::InvalidArgumentExceptionThrows if numberOfValues == 0.
OpenMPCD::NULLPointerExceptionThrows if output is nullptr.
OpenMPCD::InvalidArgumentExceptionIf OPENMPCD_DEBUG is defined, throws if values is not a Device pointer.
OpenMPCD::InvalidArgumentExceptionIf OPENMPCD_DEBUG is defined, throws if output is not a Device pointer.
Template Parameters
TThe type of values.
Parameters
[in]valuesThe values to take the average of, stored on the CUDA Device. In total, values must hold at least numberOfValues * sizeof(T) bytes.
[in]numberOfValuesThe number of values to treat.
[out]outputThe CUDA Device pointer to save the result to.

Definition at line 33 of file ImplementationDetails/Average.hpp.

◆ arithmeticMean_kernel()

template<typename T >
__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.

Template Parameters
TThe type of values.
Parameters
[in]valuesThe values to take the average of, stored on the CUDA Device. In total, values must hold at least numberOfValues * sizeof(T) bytes.
[in]numberOfValuesThe number of values to treat.
[out]outputThe CUDA Device pointer to save the result to.

Definition at line 22 of file ImplementationDetails/Average.hpp.

◆ atomicAdd() [1/2]

OPENMPCD_CUDA_DEVICE void OpenMPCD::CUDA::DeviceCode::atomicAdd ( double *const  target,
const double  increment 
)

Atomically adds increment to target.

Parameters
[in]targetThe address of the value to increment.
[in]incrementThe value to add to target.

◆ atomicAdd() [2/2]

OPENMPCD_CUDA_DEVICE void OpenMPCD::CUDA::DeviceCode::atomicAdd ( float *const  target,
const float  increment 
)

Atomically adds increment to target.

Parameters
[in]targetThe address of the value to increment.
[in]incrementThe value to add to target.

◆ collisionCellContributions()

__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.

Parameters
[in]workUnitOffsetThe number of particles to skip.
[in]particleCountThe number of particles there are.
[in]velocitiesThe particle velocities.
[in]collisionCellIndicesThe collision cell indices for the particles.
[in,out]collisionCellMomentaThe 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]collisionCellMassesThe 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]particleMassThe mass of any one particle.

◆ collisionCellStochasticRotationStep1()

__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.

Parameters
[in]workUnitOffsetThe number of particles to skip.
[in]particleCountThe number of particles there are.
[in,out]velocitiesThe 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]particleMassThe mass of any one particle.
[in]collisionCellIndicesAn array which holds, for each particle, the index of the collision cell it is currently in.
[in]collisionCellMomentaAn array holding, for each collision cell, its three cartesian momentum coordinates.
[in]collisionCellMassesAn array holding, for each collision cell, its total mass content.
[in]collisionCellRotationAxesAn array holding, for each collision cell, the random rotation axis it is to use during the collision.
[out]collisionCellFrameInternalKineticEnergiesAn 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]collisionCellParticleCountsAn 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.

◆ collisionCellStochasticRotationStep2()

__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.

Parameters
[in]workUnitOffsetThe number of particles to skip.
[in]particleCountThe number of particles there are.
[in,out]velocitiesThe 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]collisionCellIndicesAn array which holds, for each particle, the index of the collision cell it is currently in.
[in]collisionCellMomentaAn array holding, for each collision cell, its three Cartesian momentum coordinates.
[in]collisionCellMassesAn array holding, for each collision cell, its total mass content.
[in]collisionCellVelocityScalingsFor each collision cell, stores the factor by which to scale the relative velocities.

◆ constructGPURNGs()

__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.

Parameters
[in]countThe number of instances to construct.
[out]locationThe location in memory to store the instances in; memory must have been allocated at this location prior to calling this function.
[in]seedThe seed to use for the RNGs.

◆ destroyGPURNGs()

__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.

Parameters
[in]countThe number of instances to construct.
[out]locationThe location in memory where there are instances; memory must be freed manually after this function call.

◆ generateCollisionCellMBSFactors()

__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.

Parameters
[in]workUnitOffsetThe number of particles to skip.
[in]collisionCellCountThe number of collision cells there are.
[in]collisionCellFrameInternalKineticEnergiesStores, 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]collisionCellParticleCountsThe number of particles in the collision cells.
[out]collisionCellRelativeVelocityScalingsFor each collision cell, stores the factor by which to scale the relative velocities.
[in]bulkThermostatTargetkTThe target temperature for the thermostat.
[in,out]rngsThe random number generators, of which there must be at least collisionCellCount.

◆ generateCollisionCellRotationAxes()

__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.

Parameters
[in]workUnitOffsetThe number of particles to skip.
[in]collisionCellCountThe number of collision cells there are.
[out]collisionCellRotationAxesFor each collision cell, stores the Cartesian coordinates of the unit-length axis of rotation.
[in,out]rngsThe random number generators, of which there must be at least collisionCellCount.

◆ generateGridShiftVector()

__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.

Parameters
[out]gridShiftWhere to store the three grid shift coordinates to.
[in]gridShiftScaleThe scale for the grid shift vector coordinates.
[in,out]rngsPointer to at least one random number generator.

◆ getCollisionCellCenterOfMassVelocity()

__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.

Parameters
[in]collisionCellIndexThe index of the collision cell.
[in]collisionCellMomentaThe momenta of the collision cells.
[in]collisionCellMassesThe masses of the collision cells.

Definition at line 113 of file Simulation.cpp.

◆ getCollisionCellIndex()

__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, whereboxSizeX andboxSizeYare the number of collision cells in the primary simulation volume along thexandy` directions, respectively.

This function requires that OpenMPCD::CUDA::DeviceCode::setSimulationBoxSizeSymbols has been called before.

Parameters
[in]positionThe 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.

◆ getImageUnderLeesEdwardsBoundaryConditions()

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.

See also
LeesEdwardsBoundaryConditions
Parameters
[in]mpcTimeThe simulation time for the MPC steps.
[in]positionThe particle position.
[out]velocityCorrectionSet to the velocity along the x direction that needs to be added to the particle's velocity.

◆ isInPrimarySimulationVolume()

__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.

Parameters
[in]positionThe position to check.

◆ pow() [1/2]

template<typename B >
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.

Template Parameters
BThe type of the base argument, which must be an integral type.
Parameters
[in]baseThe base.
[in]exponentThe exponent.

Definition at line 50 of file Utilities.hpp.

◆ pow() [2/2]

OPENMPCD_CUDA_HOST_AND_DEVICE double OpenMPCD::CUDA::DeviceCode::pow ( const double  base,
const double  exponent 
)
inline

The power function.

Parameters
[in]baseThe base.
[in]exponentThe exponent.

Definition at line 62 of file Utilities.hpp.

◆ resetCollisionCellData()

__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.

Parameters
[in]workUnitOffsetThe number of collision cells to skip.
[in]collisionCellCountThe number of collision cells there are.
[out]collisionCellParticleCountsPointer to the memory where the collision cell particle counts are stored.
[out]collisionCellMomentaPointer to the memory where the collision cell momenta are stored (three per collision cell).
[out]collisionCellFrameInternalKineticEnergiesPointer to the memory where the collision cell kinetic energies (in the internal frame) are stored.
[out]collisionCellMassesPointer to the memory where the collision cell masses are stored.

◆ setLeesEdwardsSymbols()

void OpenMPCD::CUDA::DeviceCode::setLeesEdwardsSymbols ( const FP  shearRate,
const unsigned int  simBoxY 
)

Sets CUDA symbols for Lees-Edwards boundary conditions.

See also
LeesEdwardsBoundaryConditions

This function stores \( \Delta v_x = \dot{\gamma} L_y \) in leesEdwardsRelativeLayerVelocity.

Parameters
[in]shearRateThe shear rate \( \dot{\gamma} \).
[in]simBoxYThe size \( L_y \) of the primary simulation box along the y direction.

◆ setMPCStreamingTimestep()

void OpenMPCD::CUDA::DeviceCode::setMPCStreamingTimestep ( const FP  timestep)

Sets the symbol OpenMPCD::CUDA::DeviceCode::streamingTimestep.

Exceptions
OpenMPCD::InvalidArgumentExceptionIf OPENMPCD_DEBUG is defined, throws if timestep == 0.
Parameters
[in]timestepThe streaming timestep for each MPC fluid streaming step.

◆ setSimulationBoxSizeSymbols()

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.

Exceptions
OpenMPCD::InvalidArgumentExceptionIf OPENMPCD_DEBUG is defined, throws if either argument is 0.
Parameters
[in]simBoxSizeXThe value for mpcSimulationBoxSizeX.
[in]simBoxSizeYThe value for mpcSimulationBoxSizeY.
[in]simBoxSizeZThe value for mpcSimulationBoxSizeZ.

◆ setSRDCollisionAngleSymbol()

void OpenMPCD::CUDA::DeviceCode::setSRDCollisionAngleSymbol ( const FP  angle)

Sets the symbol OpenMPCD::CUDA::DeviceCode::srdCollisionAngle.

Parameters
[in]angleThe SRD collision angle.

◆ sortIntoCollisionCellsLeesEdwards()

__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.

See also
LeesEdwardsBoundaryConditions
Parameters
[in]particleIDThe particle ID.
[in]gridShift_The grid shift vector to temporarily add to the position.
[in]mpcTimeThe MPC time that has passed since the start of the simulation.
[in]positionsThe particle positions.
[in,out]velocitiesThe particle velocities, which may be changed due to Lees-Edwards boundary conditions.
[out]velocityCorrectionsThe velocity corrections applied due to Lees-Edwards boundary conditions.
[out]collisionCellIndicesThe collision cell indices for the particles.

◆ sortParticlesIntoCollisionCellsLeesEdwards()

__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.

See also
LeesEdwardsBoundaryConditions
Parameters
[in]workUnitOffsetThe number of particles to skip.
[in]particleCountThe number of particles there are.
[in]gridShift_The grid shift vector to temporarily add to all positions.
[in]mpcTimeThe MPC time that has passed since the start of the simulation.
[in]positionsThe MPC fluid particle positions.
[in,out]velocitiesThe MPC fluid particle velocities, which may be changed due to Lees-Edwards boundary conditions.
[out]velocityCorrectionsThe 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]collisionCellIndicesThe collision cell indices for the particles.

◆ undoLeesEdwardsVelocityCorrections()

__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.

Parameters
[in]workUnitOffsetThe number of particles to skip.
[in]particleCountThe number of particles there are.
[in,out]velocitiesThe particle velocities.
[in]velocityCorrectionsThe velocity corrections as returned by sortIntoCollisionCellsLeesEdwards.

◆ velocityVerletStep1() [1/2]

__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.

Parameters
[in]positionThe position.
[in]velocityThe velocity.
[in]accelerationThe acceleration.
[in]timestepThe timestep.

◆ velocityVerletStep1() [2/2]

__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.

Parameters
[in,out]positionThe position.
[in]velocityThe velocity.
[in]accelerationThe acceleration.
[in]timestepThe timestep.

◆ velocityVerletStep2() [1/2]

__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.

Parameters
[in]velocityThe velocity.
[in]oldAccelerationThe acceleration at the time \(t\).
[in]newAccelerationThe acceleration at the time \(t+h\).
[in]timestepThe timestep.

◆ velocityVerletStep2() [2/2]

__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.

Parameters
[in,out]velocityThe velocity.
[in]oldAccelerationThe acceleration at the time \(t\).
[in]newAccelerationThe acceleration at the time \(t+h\).
[in]timestepThe timestep.

Variable Documentation

◆ collisionCellCount

__constant__ unsigned int OpenMPCD::CUDA::DeviceCode::collisionCellCount

The number of collision cells in the system.

◆ leesEdwardsRelativeLayerVelocity

__constant__ FP OpenMPCD::CUDA::DeviceCode::leesEdwardsRelativeLayerVelocity

The relative x-velocity of y-adjacent layers in Lees-Edwards boundary conditions.

See also
LeesEdwardsBoundaryConditions

◆ mpcSimulationBoxSizeX

__constant__ unsigned int OpenMPCD::CUDA::DeviceCode::mpcSimulationBoxSizeX

The size of the primary simulation box along the x direction.

◆ mpcSimulationBoxSizeY

__constant__ unsigned int OpenMPCD::CUDA::DeviceCode::mpcSimulationBoxSizeY

The size of the primary simulation box along the y direction.

◆ mpcSimulationBoxSizeZ

__constant__ unsigned int OpenMPCD::CUDA::DeviceCode::mpcSimulationBoxSizeZ

The size of the primary simulation box along the z direction.

◆ srdCollisionAngle

__constant__ FP OpenMPCD::CUDA::DeviceCode::srdCollisionAngle

The collision angle for SRD collisions.

◆ streamingTimestep

__constant__ FP OpenMPCD::CUDA::DeviceCode::streamingTimestep

The MPC streaming time-step.