OpenMPCD
DeviceCode/Chains.hpp
Go to the documentation of this file.
1 /**
2  * @file
3  * Defines CUDA Device code for MPC fluids consisting out of MPC-particle chains.
4  */
5 
6 #ifndef OPENMPCD_CUDA_MPCFLUID_DEVICECODE_CHAINS_HPP
7 #define OPENMPCD_CUDA_MPCFLUID_DEVICECODE_CHAINS_HPP
8 
9 #include <OpenMPCD/Types.hpp>
10 #include <OpenMPCD/Vector3D.hpp>
11 
12 namespace OpenMPCD
13 {
14 namespace CUDA
15 {
16 namespace MPCFluid
17 {
18 namespace DeviceCode
19 {
20  /**
21  * Saves the center-of-mass velocities of the MPC fluid's chains to the given buffer.
22  *
23  * This function assumes that all chains have the same number of particles,
24  * and that all particles have the same mass.
25  *
26  * This function requires that the
27  * `OpenMPCD::CUDA::MPCFluid::DeviceCode::mpcParticleCount`
28  * symbol has been set, and its value is equal to the `mpcParticleCount`
29  * passed to this function.
30  *
31  * This function assumes that `mpcParticleCount` is an integer multiple of
32  * `chainLength`, and assumes both of them to be non-zero.
33  *
34  * @param[in] mpcParticleCount The number of individual MPC fluid particles.
35  * @param[in] chainLength The number of individual MPC fluid particles per chain.
36  * @param[in] velocities The Device array of MPC fluid particle velocities.
37  * @param[out] comVelocities Device buffer where the center-of-mass velocities are saved.
38  */
40  const unsigned int mpcParticleCount,
41  const unsigned int chainLength,
42  const MPCParticleVelocityType* const velocities,
43  MPCParticleVelocityType* const comVelocities);
44 
45  /**
46  * Saves the center-of-mass velocities of the MPC fluid's chains to the given buffer.
47  *
48  * This function assumes that all chains have the same number of particles,
49  * and that all particles have the same mass.
50  *
51  * This function requires that the
52  * `OpenMPCD::CUDA::MPCFluid::DeviceCode::mpcParticleCount`
53  * symbol has been set to a value that is non-zero and an integer multiple
54  * of `chainLength`.
55  *
56  * @param[in] workUnitOffset The number of chains to skip.
57  * @param[in] chainLength The number of individual MPC fluid particles per chain.
58  * @param[in] velocities The array of MPC fluid particle velocities.
59  * @param[out] comVelocities Buffer where the center-of-mass velocities are saved.
60  */
62  const unsigned int workUnitOffset,
63  const unsigned int chainLength,
64  const MPCParticleVelocityType* const velocities,
65  MPCParticleVelocityType* const comVelocities);
66 
67  /**
68  * Returns the center-of-mass velocity of for the given chain.
69  *
70  * This function assumes that all chains have the same number of particles,
71  * and that all particles have the same mass.
72  *
73  * @param[in] chainID The ID of the chain.
74  * @param[in] chainLength The number of individual MPC fluid particles per chain.
75  * @param[in] velocities The array of MPC fluid particle velocities.
76  */
77  __device__ Vector3D<MPCParticleVelocityType> getCenterOfMassVelocity_chain(
78  const unsigned int chainID,
79  const unsigned int chainLength,
80  const MPCParticleVelocityType* const velocities);
81 } //namespace DeviceCode
82 } //namespace MPCFluid
83 } //namespace CUDA
84 } //namespace OpenMPCD
85 
86 #endif
OpenMPCD::CUDA::MPCFluid::DeviceCode::getCenterOfMassVelocities_chain
void getCenterOfMassVelocities_chain(const unsigned int mpcParticleCount, const unsigned int chainLength, const MPCParticleVelocityType *const velocities, MPCParticleVelocityType *const comVelocities)
Saves the center-of-mass velocities of the MPC fluid's chains to the given buffer.
OpenMPCD::CUDA::MPCFluid::DeviceCode::getCenterOfMassVelocities_chain_kernel
__global__ void getCenterOfMassVelocities_chain_kernel(const unsigned int workUnitOffset, const unsigned int chainLength, const MPCParticleVelocityType *const velocities, MPCParticleVelocityType *const comVelocities)
Saves the center-of-mass velocities of the MPC fluid's chains to the given buffer.
OpenMPCD::CUDA::MPCFluid::DeviceCode::mpcParticleCount
__constant__ unsigned int mpcParticleCount
The number of MPC fluid particles.
OpenMPCD::CUDA::MPCFluid::DeviceCode::getCenterOfMassVelocity_chain
__device__ Vector3D< MPCParticleVelocityType > getCenterOfMassVelocity_chain(const unsigned int chainID, const unsigned int chainLength, const MPCParticleVelocityType *const velocities)
Returns the center-of-mass velocity of for the given chain.
Vector3D.hpp
Types.hpp