Fix random typos
[alexxy/gromacs.git] / src / gromacs / gpu_utils / vectype_ops.cuh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 #ifndef VECTYPE_OPS_CUH
37 #define VECTYPE_OPS_CUH
38
39 /**** float3 ****/
40 __forceinline__ __host__ __device__ float3 make_float3(float s)
41 {
42     return make_float3(s, s, s);
43 }
44 __forceinline__ __host__ __device__ float3 make_float3(float4 a)
45 {
46     return make_float3(a.x, a.y, a.z);
47 }
48 __forceinline__ __host__ __device__ float3 operator-(const float3& a)
49 {
50     return make_float3(-a.x, -a.y, -a.z);
51 }
52 __forceinline__ __host__ __device__ float3 operator+(float3 a, float3 b)
53 {
54     return make_float3(a.x + b.x, a.y + b.y, a.z + b.z);
55 }
56 __forceinline__ __host__ __device__ float3 operator-(float3 a, float3 b)
57 {
58     return make_float3(a.x - b.x, a.y - b.y, a.z - b.z);
59 }
60 __forceinline__ __host__ __device__ float3 operator*(float3 a, float k)
61 {
62     return make_float3(k * a.x, k * a.y, k * a.z);
63 }
64 __forceinline__ __host__ __device__ float3 operator*(float k, float3 a)
65 {
66     return make_float3(k * a.x, k * a.y, k * a.z);
67 }
68 // NOLINTNEXTLINE(google-runtime-references)
69 __forceinline__ __host__ __device__ void operator+=(float3& a, float3 b)
70 {
71     a.x += b.x;
72     a.y += b.y;
73     a.z += b.z;
74 }
75 // NOLINTNEXTLINE(google-runtime-references)
76 __forceinline__ __host__ __device__ void operator+=(float3& a, float4 b)
77 {
78     a.x += b.x;
79     a.y += b.y;
80     a.z += b.z;
81 }
82 // NOLINTNEXTLINE(google-runtime-references)
83 __forceinline__ __host__ __device__ void operator-=(float3& a, float3 b)
84 {
85     a.x -= b.x;
86     a.y -= b.y;
87     a.z -= b.z;
88 }
89 __forceinline__ __host__ __device__ float norm(float3 a)
90 {
91     return sqrt(a.x * a.x + a.y * a.y + a.z * a.z);
92 }
93 __forceinline__ __host__ __device__ float norm2(float3 a)
94 {
95     return (a.x * a.x + a.y * a.y + a.z * a.z);
96 }
97 __forceinline__ __host__ __device__ float dist3(float3 a, float3 b)
98 {
99     return norm(b - a);
100 }
101 __forceinline__ __host__ __device__ float3 operator*(float3 a, float3 b)
102 {
103     return make_float3(a.x * b.x, a.y * b.y, a.z * b.z);
104 }
105 // NOLINTNEXTLINE(google-runtime-references)
106 __forceinline__ __host__ __device__ void operator*=(float3& a, float3 b)
107 {
108     a.x *= b.x;
109     a.y *= b.y;
110     a.z *= b.z;
111 }
112 // NOLINTNEXTLINE(google-runtime-references)
113 __forceinline__ __host__ __device__ void operator*=(float3& a, float b)
114 {
115     a.x *= b;
116     a.y *= b;
117     a.z *= b;
118 }
119 __forceinline__ __device__ void atomicAdd(float3* addr, float3 val)
120 {
121     atomicAdd(&addr->x, val.x);
122     atomicAdd(&addr->y, val.y);
123     atomicAdd(&addr->z, val.z);
124 }
125 /****************************************************************/
126
127 /**** float4 ****/
128 __forceinline__ __host__ __device__ float4 make_float4(float s)
129 {
130     return make_float4(s, s, s, s);
131 }
132 __forceinline__ __host__ __device__ float4 make_float4(float3 a)
133 {
134     return make_float4(a.x, a.y, a.z, 0.0F);
135 }
136 __forceinline__ __host__ __device__ float4 operator+(float4 a, float4 b)
137 {
138     return make_float4(a.x + b.x, a.y + b.y, a.z + b.z, a.w + b.w);
139 }
140 __forceinline__ __host__ __device__ float4 operator+(float4 a, float3 b)
141 {
142     return make_float4(a.x + b.x, a.y + b.y, a.z + b.z, a.w);
143 }
144 __forceinline__ __host__ __device__ float4 operator-(float4 a, float4 b)
145 {
146     return make_float4(a.x - b.x, a.y - b.y, a.z - b.z, a.w - b.w);
147 }
148 __forceinline__ __host__ __device__ float4 operator*(float4 a, float k)
149 {
150     return make_float4(k * a.x, k * a.y, k * a.z, k * a.w);
151 }
152 __forceinline__ __host__ __device__ void operator+=(float4& a, float4 b)
153 {
154     a.x += b.x;
155     a.y += b.y;
156     a.z += b.z;
157     a.w += b.w;
158 }
159 // NOLINTNEXTLINE(google-runtime-references)
160 __forceinline__ __host__ __device__ void operator+=(float4& a, float3 b)
161 {
162     a.x += b.x;
163     a.y += b.y;
164     a.z += b.z;
165 }
166 // NOLINTNEXTLINE(google-runtime-references)
167 __forceinline__ __host__ __device__ void operator-=(float4& a, float3 b)
168 {
169     a.x -= b.x;
170     a.y -= b.y;
171     a.z -= b.z;
172 }
173
174 __forceinline__ __host__ __device__ float norm(float4 a)
175 {
176     return sqrt(a.x * a.x + a.y * a.y + a.z * a.z + a.w * a.w);
177 }
178
179 __forceinline__ __host__ __device__ float dist3(float4 a, float4 b)
180 {
181     return norm(b - a);
182 }
183
184 /* \brief Compute the scalar product of two vectors.
185  *
186  * \param[in] a  First vector.
187  * \param[in] b  Second vector.
188  * \returns Scalar product.
189  */
190 __forceinline__ __device__ float iprod(const float3 a, const float3 b)
191 {
192     return a.x * b.x + a.y * b.y + a.z * b.z;
193 }
194
195 /* \brief Compute the vector product of two vectors.
196  *
197  * \param[in] a  First vector.
198  * \param[in] b  Second vector.
199  * \returns Vector product.
200  */
201 __forceinline__ __device__ float3 cprod(const float3 a, const float3 b)
202 {
203     float3 c;
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;
207     return c;
208 }
209
210 /* \brief Cosine of an angle between two vectors.
211  *
212  * Computes cosine using the following formula:
213  *
214  *                  ax*bx + ay*by + az*bz
215  * cos-vec (a,b) =  ---------------------
216  *                      ||a|| * ||b||
217  *
218  * This function also makes sure that the cosine does not leave the [-1, 1]
219  * interval, which can happen due to numerical errors.
220  *
221  * \param[in] a  First vector.
222  * \param[in] b  Second vector.
223  * \returns Cosine between a and b.
224  */
225 __forceinline__ __device__ float cos_angle(const float3 a, const float3 b)
226 {
227     float cosval;
228
229     float ipa  = norm2(a);
230     float ipb  = norm2(b);
231     float ip   = iprod(a, b);
232     float ipab = ipa * ipb;
233     if (ipab > 0.0F)
234     {
235         cosval = ip * rsqrt(ipab);
236     }
237     else
238     {
239         cosval = 1.0F;
240     }
241     if (cosval > 1.0F)
242     {
243         return 1.0F;
244     }
245     if (cosval < -1.0F)
246     {
247         return -1.0F;
248     }
249
250     return cosval;
251 }
252
253 /* \brief Compute the angle between two vectors.
254  *
255  * Uses atan( |axb| / a.b ) formula.
256  *
257  * \param[in] a  First vector.
258  * \param[in] b  Second vector.
259  * \returns Angle between vectors in radians.
260  */
261 __forceinline__ __device__ float gmx_angle(const float3 a, const float3 b)
262 {
263     float3 w = cprod(a, b);
264
265     float wlen = norm(w);
266     float s    = iprod(a, b);
267
268     return atan2f(wlen, s); // requires float
269 }
270
271 /* \brief Atomically add components of the vector.
272  *
273  * Executes atomicAdd one-by-one on all components of the float3 vector.
274  *
275  * \param[in] a  First vector.
276  * \param[in] b  Second vector.
277  */
278 // NOLINTNEXTLINE(google-runtime-references)
279 __forceinline__ __device__ void atomicAdd(float3& a, const float3 b)
280 {
281     atomicAdd(&a.x, b.x);
282     atomicAdd(&a.y, b.y);
283     atomicAdd(&a.z, b.z);
284 }
285
286 #endif /* VECTYPE_OPS_CUH */