Fix error in ensureReferenceCount
[gromacs.git] / src / gromacs / gpu_utils / vectype_ops.cuh
blobcce3fc90082e00a5e4730c9d980a746858f2ea5b
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2015,2016,2019,2020, 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  */
36 #ifndef VECTYPE_OPS_CUH
37 #define VECTYPE_OPS_CUH
39 /**** float3 ****/
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-(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 __forceinline__ __host__ __device__ void operator+=(float3& a, float3 b)
70     a.x += b.x;
71     a.y += b.y;
72     a.z += b.z;
74 __forceinline__ __host__ __device__ void operator+=(float3& a, float4 b)
76     a.x += b.x;
77     a.y += b.y;
78     a.z += b.z;
80 __forceinline__ __host__ __device__ void operator-=(float3& a, float3 b)
82     a.x -= b.x;
83     a.y -= b.y;
84     a.z -= b.z;
86 __forceinline__ __host__ __device__ float norm(float3 a)
88     return sqrt(a.x * a.x + a.y * a.y + a.z * a.z);
90 __forceinline__ __host__ __device__ float norm2(float3 a)
92     return (a.x * a.x + a.y * a.y + a.z * a.z);
94 __forceinline__ __host__ __device__ float dist3(float3 a, float3 b)
96     return norm(b - a);
98 __forceinline__ __host__ __device__ float3 operator*(float3 a, float3 b)
100     return make_float3(a.x * b.x, a.y * b.y, a.z * b.z);
102 __forceinline__ __host__ __device__ void operator*=(float3& a, float3 b)
104     a.x *= b.x;
105     a.y *= b.y;
106     a.z *= b.z;
108 __forceinline__ __host__ __device__ void operator*=(float3& a, float b)
110     a.x *= b;
111     a.y *= b;
112     a.z *= b;
114 __forceinline__ __device__ void atomicAdd(float3* addr, float3 val)
116     atomicAdd(&addr->x, val.x);
117     atomicAdd(&addr->y, val.y);
118     atomicAdd(&addr->z, val.z);
120 /****************************************************************/
122 /**** float4 ****/
123 __forceinline__ __host__ __device__ float4 make_float4(float s)
125     return make_float4(s, s, s, s);
127 __forceinline__ __host__ __device__ float4 make_float4(float3 a)
129     return make_float4(a.x, a.y, a.z, 0.0f);
131 __forceinline__ __host__ __device__ float4 operator+(float4 a, float4 b)
133     return make_float4(a.x + b.x, a.y + b.y, a.z + b.z, a.w + b.w);
135 __forceinline__ __host__ __device__ float4 operator+(float4 a, float3 b)
137     return make_float4(a.x + b.x, a.y + b.y, a.z + b.z, a.w);
139 __forceinline__ __host__ __device__ float4 operator-(float4 a, float4 b)
141     return make_float4(a.x - b.x, a.y - b.y, a.z - b.z, a.w - b.w);
143 __forceinline__ __host__ __device__ float4 operator*(float4 a, float k)
145     return make_float4(k * a.x, k * a.y, k * a.z, k * a.w);
147 __forceinline__ __host__ __device__ void operator+=(float4& a, float4 b)
149     a.x += b.x;
150     a.y += b.y;
151     a.z += b.z;
152     a.w += b.w;
154 __forceinline__ __host__ __device__ void operator+=(float4& a, float3 b)
156     a.x += b.x;
157     a.y += b.y;
158     a.z += b.z;
160 __forceinline__ __host__ __device__ void operator-=(float4& a, float3 b)
162     a.x -= b.x;
163     a.y -= b.y;
164     a.z -= b.z;
167 __forceinline__ __host__ __device__ float norm(float4 a)
169     return sqrt(a.x * a.x + a.y * a.y + a.z * a.z + a.w * a.w);
172 __forceinline__ __host__ __device__ float dist3(float4 a, float4 b)
174     return norm(b - a);
177 /* \brief Compute the scalar product of two vectors.
179  * \param[in] a  First vector.
180  * \param[in] b  Second vector.
181  * \returns Scalar product.
182  */
183 __forceinline__ __device__ float iprod(const float3 a, const float3 b)
185     return a.x * b.x + a.y * b.y + a.z * b.z;
188 /* \brief Compute the vector product of two vectors.
190  * \param[in] a  First vector.
191  * \param[in] b  Second vector.
192  * \returns Vector product.
193  */
194 __forceinline__ __device__ float3 cprod(const float3 a, const float3 b)
196     float3 c;
197     c.x = a.y * b.z - a.z * b.y;
198     c.y = a.z * b.x - a.x * b.z;
199     c.z = a.x * b.y - a.y * b.x;
200     return c;
203 /* \brief Cosine of an angle between two vectors.
205  * Computes cosine using the following formula:
207  *                  ax*bx + ay*by + az*bz
208  * cos-vec (a,b) =  ---------------------
209  *                      ||a|| * ||b||
211  * This function also makes sure that the cosine does not leave the [-1, 1]
212  * interval, which can happen due to numerical errors.
214  * \param[in] a  First vector.
215  * \param[in] b  Second vector.
216  * \returns Cosine between a and b.
217  */
218 __forceinline__ __device__ float cos_angle(const float3 a, const float3 b)
220     float cosval;
222     float ipa  = norm2(a);
223     float ipb  = norm2(b);
224     float ip   = iprod(a, b);
225     float ipab = ipa * ipb;
226     if (ipab > 0.0f)
227     {
228         cosval = ip * rsqrt(ipab);
229     }
230     else
231     {
232         cosval = 1.0f;
233     }
234     if (cosval > 1.0f)
235     {
236         return 1.0f;
237     }
238     if (cosval < -1.0f)
239     {
240         return -1.0f;
241     }
243     return cosval;
246 /* \brief Compute the angle between two vectors.
248  * Uses atan( |axb| / a.b ) formula.
250  * \param[in] a  First vector.
251  * \param[in] b  Second vector.
252  * \returns Angle between vectors in radians.
253  */
254 __forceinline__ __device__ float gmx_angle(const float3 a, const float3 b)
256     float3 w = cprod(a, b);
258     float wlen = norm(w);
259     float s    = iprod(a, b);
261     return atan2f(wlen, s); // requires float
264 /* \brief Atomically add components of the vector.
266  * Executes atomicAdd one-by-one on all components of the float3 vector.
268  * \param[in] a  First vector.
269  * \param[in] b  Second vector.
270  * \returns Angle between vectors.
271  */
272 __forceinline__ __device__ void atomicAdd(float3& a, const float3 b)
274     atomicAdd(&a.x, b.x);
275     atomicAdd(&a.y, b.y);
276     atomicAdd(&a.z, b.z);
279 #endif /* VECTYPE_OPS_CUH */