OpenMPCD
CUDA/DeviceCode/ImplementationDetails/NormalMode.hpp
Go to the documentation of this file.
1 /**
2  * @file
3  * Implements functionality in the `OpenMPCD::CUDA::DeviceCode::NormalMode`
4  * namespace.
5  */
6 
7 #ifndef OPENMPCD_CUDA_DEVICECODE_IMPLEMENTATIONDETAILS_NORMALMODE_HPP
8 #define OPENMPCD_CUDA_DEVICECODE_IMPLEMENTATIONDETAILS_NORMALMODE_HPP
9 
13 
14 #include <boost/static_assert.hpp>
15 #include <boost/type_traits/is_floating_point.hpp>
16 
17 namespace OpenMPCD
18 {
19 namespace CUDA
20 {
21 namespace DeviceCode
22 {
23 namespace NormalMode
24 {
25 
26 template<typename T>
29  const unsigned int i, const Vector3D<T>* const vectors,
30  const std::size_t N, const T shift)
31 {
32  BOOST_STATIC_ASSERT(boost::is_floating_point<T>::value);
33  //non-floating-point `T` is probably a mistake
34 
35  OPENMPCD_DEBUG_ASSERT(vectors);
36  OPENMPCD_DEBUG_ASSERT(N != 0);
37  OPENMPCD_DEBUG_ASSERT(i <= N);
38 
39  Vector3D<T> ret(0, 0, 0);
40 
41  const T argPart = T(i) / T(N);
42  for(std::size_t n = 1; n <= N; ++n)
43  {
44  const T arg = argPart * (n + shift);
45  ret += vectors[n - 1] * Utility::MathematicalFunctions::cospi(arg);
46  }
47 
48  return ret / T(N);
49 }
50 
51 template<typename T>
54  const unsigned int i, const T* const vectors, const std::size_t N,
55  const T shift)
56 {
57  BOOST_STATIC_ASSERT(boost::is_floating_point<T>::value);
58  //non-floating-point `T` is probably a mistake
59 
60  OPENMPCD_DEBUG_ASSERT(vectors);
61  OPENMPCD_DEBUG_ASSERT(N != 0);
62  OPENMPCD_DEBUG_ASSERT(i <= N);
63 
64  Vector3D<T> ret(0, 0, 0);
65 
66  const T argPart = T(i) / T(N);
67  for(std::size_t n = 1; n <= N; ++n)
68  {
69  const OpenMPCD::RemotelyStoredVector<const T> v(vectors, n - 1);
70  ret +=
71  v * Utility::MathematicalFunctions::cospi(argPart * (n + shift));
72  }
73 
74  return ret / T(N);
75 }
76 
77 template<typename T>
80  const T* const vectors, const std::size_t N, T* const result,
81  const T shift)
82 {
83  BOOST_STATIC_ASSERT(boost::is_floating_point<T>::value);
84  //non-floating-point `T` is probably a mistake
85 
86  for(std::size_t i = 0; i <= N; ++i)
87  {
89  r = computeNormalCoordinate(i, vectors, N, shift);
90  }
91 }
92 
93 template<typename T>
94 __global__ void computeNormalCoordinates(
95  const unsigned int workUnitOffset,
96  const unsigned int chainLength,
97  const unsigned int chainCount,
98  const T* const positions,
99  T* const normalModeCoordinates,
100  const T shift)
101 {
102  BOOST_STATIC_ASSERT(boost::is_floating_point<T>::value);
103  //non-floating-point `T` is probably a mistake
104 
105  OPENMPCD_DEBUG_ASSERT(chainLength != 0);
106  OPENMPCD_DEBUG_ASSERT(positions);
107  OPENMPCD_DEBUG_ASSERT(normalModeCoordinates);
108 
109  const unsigned int chainID =
110  blockIdx.x * blockDim.x + threadIdx.x + workUnitOffset;
111 
112  if(chainID >= chainCount)
113  return;
114 
116  positions + 3 * chainID * chainLength,
117  chainLength,
118  normalModeCoordinates + 3 * chainID * (chainLength + 1),
119  shift);
120 }
121 
122 } //namespace NormalMode
123 } //namespace DeviceCode
124 } //namespace CUDA
125 } //namespace OpenMPCD
126 
127 #endif //OPENMPCD_CUDA_DEVICECODE_IMPLEMENTATIONDETAILS_NORMALMODE_HPP
OpenMPCD::CUDA::DeviceCode::NormalMode::computeNormalCoordinates
OPENMPCD_CUDA_DEVICE void computeNormalCoordinates(const T *const vectors, const std::size_t N, T *const result, const T shift)
Calculates all normal coordinates for the given array of vectors.
Definition: CUDA/DeviceCode/ImplementationDetails/NormalMode.hpp:79
OpenMPCD::RemotelyStoredVector
Represents a vector whose data is stored elsewhere.
Definition: RemotelyStoredVector.hpp:26
OpenMPCD::Vector3D
3-dimensional vector.
Definition: Vector3D.hpp:38
RemotelyStoredVector.hpp
OPENMPCD_DEBUG_ASSERT
#define OPENMPCD_DEBUG_ASSERT(assertion)
Asserts that the given expression evaluates to true, but only if OPENMPCD_DEBUG is defined.
Definition: OPENMPCD_DEBUG_ASSERT.hpp:88
OpenMPCD::Utility::MathematicalFunctions::cospi
OPENMPCD_CUDA_HOST_AND_DEVICE T cospi(const T x)
Returns the cosine of the product of the argument and .
OPENMPCD_CUDA_DEVICE
#define OPENMPCD_CUDA_DEVICE
Denotes a function to be callable from a CUDA Device.
Definition: Macros.hpp:25
OPENMPCD_DEBUG_ASSERT.hpp
OpenMPCD::CUDA::DeviceCode::NormalMode::computeNormalCoordinate
const OPENMPCD_CUDA_DEVICE Vector3D< T > computeNormalCoordinate(const unsigned int i, const Vector3D< T > *const vectors, const std::size_t N, const T shift)
Calculates a normal coordinate.
Definition: CUDA/DeviceCode/ImplementationDetails/NormalMode.hpp:28
MathematicalFunctions.hpp