2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2015,2016,2019,2020,2021, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 #ifndef VECTYPE_OPS_CUH
37 #define VECTYPE_OPS_CUH
40 __forceinline__ __host__ __device__ float3 make_float3(float s)
42 return make_float3(s, s, s);
44 __forceinline__ __host__ __device__ float3 make_float3(float4 a)
46 return make_float3(a.x, a.y, a.z);
48 __forceinline__ __host__ __device__ float3 operator-(const float3& a)
50 return make_float3(-a.x, -a.y, -a.z);
52 __forceinline__ __host__ __device__ float3 operator+(float3 a, float3 b)
54 return make_float3(a.x + b.x, a.y + b.y, a.z + b.z);
56 __forceinline__ __host__ __device__ float3 operator-(float3 a, float3 b)
58 return make_float3(a.x - b.x, a.y - b.y, a.z - b.z);
60 __forceinline__ __host__ __device__ float3 operator*(float3 a, float k)
62 return make_float3(k * a.x, k * a.y, k * a.z);
64 __forceinline__ __host__ __device__ float3 operator*(float k, float3 a)
66 return make_float3(k * a.x, k * a.y, k * a.z);
68 // NOLINTNEXTLINE(google-runtime-references)
69 __forceinline__ __host__ __device__ void operator+=(float3& a, float3 b)
75 // NOLINTNEXTLINE(google-runtime-references)
76 __forceinline__ __host__ __device__ void operator+=(float3& a, float4 b)
82 // NOLINTNEXTLINE(google-runtime-references)
83 __forceinline__ __host__ __device__ void operator-=(float3& a, float3 b)
89 __forceinline__ __host__ __device__ float norm(float3 a)
91 return sqrt(a.x * a.x + a.y * a.y + a.z * a.z);
93 __forceinline__ __host__ __device__ float norm2(float3 a)
95 return (a.x * a.x + a.y * a.y + a.z * a.z);
97 __forceinline__ __host__ __device__ float dist3(float3 a, float3 b)
101 __forceinline__ __host__ __device__ float3 operator*(float3 a, float3 b)
103 return make_float3(a.x * b.x, a.y * b.y, a.z * b.z);
105 // NOLINTNEXTLINE(google-runtime-references)
106 __forceinline__ __host__ __device__ void operator*=(float3& a, float3 b)
112 // NOLINTNEXTLINE(google-runtime-references)
113 __forceinline__ __host__ __device__ void operator*=(float3& a, float b)
119 __forceinline__ __device__ void atomicAdd(float3* addr, float3 val)
121 atomicAdd(&addr->x, val.x);
122 atomicAdd(&addr->y, val.y);
123 atomicAdd(&addr->z, val.z);
125 /****************************************************************/
128 __forceinline__ __host__ __device__ float4 make_float4(float s)
130 return make_float4(s, s, s, s);
132 __forceinline__ __host__ __device__ float4 make_float4(float3 a)
134 return make_float4(a.x, a.y, a.z, 0.0F);
136 __forceinline__ __host__ __device__ float4 operator+(float4 a, float4 b)
138 return make_float4(a.x + b.x, a.y + b.y, a.z + b.z, a.w + b.w);
140 __forceinline__ __host__ __device__ float4 operator+(float4 a, float3 b)
142 return make_float4(a.x + b.x, a.y + b.y, a.z + b.z, a.w);
144 __forceinline__ __host__ __device__ float4 operator-(float4 a, float4 b)
146 return make_float4(a.x - b.x, a.y - b.y, a.z - b.z, a.w - b.w);
148 __forceinline__ __host__ __device__ float4 operator*(float4 a, float k)
150 return make_float4(k * a.x, k * a.y, k * a.z, k * a.w);
152 __forceinline__ __host__ __device__ void operator+=(float4& a, float4 b)
159 // NOLINTNEXTLINE(google-runtime-references)
160 __forceinline__ __host__ __device__ void operator+=(float4& a, float3 b)
166 // NOLINTNEXTLINE(google-runtime-references)
167 __forceinline__ __host__ __device__ void operator-=(float4& a, float3 b)
174 __forceinline__ __host__ __device__ float norm(float4 a)
176 return sqrt(a.x * a.x + a.y * a.y + a.z * a.z + a.w * a.w);
179 __forceinline__ __host__ __device__ float dist3(float4 a, float4 b)
184 /* \brief Compute the scalar product of two vectors.
186 * \param[in] a First vector.
187 * \param[in] b Second vector.
188 * \returns Scalar product.
190 __forceinline__ __device__ float iprod(const float3 a, const float3 b)
192 return a.x * b.x + a.y * b.y + a.z * b.z;
195 /* \brief Compute the vector product of two vectors.
197 * \param[in] a First vector.
198 * \param[in] b Second vector.
199 * \returns Vector product.
201 __forceinline__ __device__ float3 cprod(const float3 a, const float3 b)
204 c.x = a.y * b.z - a.z * b.y;
205 c.y = a.z * b.x - a.x * b.z;
206 c.z = a.x * b.y - a.y * b.x;
210 /* \brief Cosine of an angle between two vectors.
212 * Computes cosine using the following formula:
214 * ax*bx + ay*by + az*bz
215 * cos-vec (a,b) = ---------------------
218 * This function also makes sure that the cosine does not leave the [-1, 1]
219 * interval, which can happen due to numerical errors.
221 * \param[in] a First vector.
222 * \param[in] b Second vector.
223 * \returns Cosine between a and b.
225 __forceinline__ __device__ float cos_angle(const float3 a, const float3 b)
229 float ipa = norm2(a);
230 float ipb = norm2(b);
231 float ip = iprod(a, b);
232 float ipab = ipa * ipb;
235 cosval = ip * rsqrt(ipab);
253 /* \brief Compute the angle between two vectors.
255 * Uses atan( |axb| / a.b ) formula.
257 * \param[in] a First vector.
258 * \param[in] b Second vector.
259 * \returns Angle between vectors in radians.
261 __forceinline__ __device__ float gmx_angle(const float3 a, const float3 b)
263 float3 w = cprod(a, b);
265 float wlen = norm(w);
266 float s = iprod(a, b);
268 return atan2f(wlen, s); // requires float
271 /* \brief Atomically add components of the vector.
273 * Executes atomicAdd one-by-one on all components of the float3 vector.
275 * \param[in] a First vector.
276 * \param[in] b Second vector.
278 // NOLINTNEXTLINE(google-runtime-references)
279 __forceinline__ __device__ void atomicAdd(float3& a, const float3 b)
281 atomicAdd(&a.x, b.x);
282 atomicAdd(&a.y, b.y);
283 atomicAdd(&a.z, b.z);
286 #endif /* VECTYPE_OPS_CUH */