2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, 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.
35 #ifndef GMX_SIMD_SIMD_MATH_H_
36 #define GMX_SIMD_SIMD_MATH_H_
38 /*! \libinternal \file
40 * \brief Math functions for SIMD datatypes.
42 * \attention This file is generic for all SIMD architectures, so you cannot
43 * assume that any of the optional SIMD features (as defined in simd.h) are
44 * present. In particular, this means you cannot assume support for integers,
45 * logical operations (neither on floating-point nor integer values), shifts,
46 * and the architecture might only have SIMD for either float or double.
47 * Second, to keep this file clean and general, any additions to this file
48 * must work for all possible SIMD architectures in both single and double
49 * precision (if they support it), and you cannot make any assumptions about
52 * \author Erik Lindahl <erik.lindahl@scilifelab.se>
55 * \ingroup module_simd
60 #include "gromacs/math/utilities.h"
61 #include "gromacs/simd/simd.h"
66 /*! \addtogroup module_simd */
69 /*! \name Implementation accuracy settings
73 /*! \brief We accept lsb errors for 1/sqrt(x) and 1/x, so float target is 22 bits */
74 #define GMX_SIMD_MATH_TARGET_SINGLE_BITS 22
76 /*! \brief We accept "double" that has 2x single precision - 44 bits.
78 * This way two Newton-Raphson iterations will suffice in double precision.
80 #define GMX_SIMD_MATH_TARGET_DOUBLE_BITS 44
84 #ifdef GMX_SIMD_HAVE_FLOAT
86 /*! \name Single precision SIMD math functions
88 * \note In most cases you should use the real-precision functions instead.
92 /****************************************
93 * SINGLE PRECISION SIMD MATH FUNCTIONS *
94 ****************************************/
96 /*! \brief SIMD float utility to sum a+b+c+d.
98 * You should normally call the real-precision routine \ref gmx_simd_sum4_r.
100 * \param a term 1 (multiple values)
101 * \param b term 2 (multiple values)
102 * \param c term 3 (multiple values)
103 * \param d term 4 (multiple values)
104 * \return sum of terms 1-4 (multiple values)
106 static gmx_inline gmx_simd_float_t gmx_simdcall
107 gmx_simd_sum4_f(gmx_simd_float_t a, gmx_simd_float_t b,
108 gmx_simd_float_t c, gmx_simd_float_t d)
110 return gmx_simd_add_f(gmx_simd_add_f(a, b), gmx_simd_add_f(c, d));
113 /*! \brief Return -a if b is negative, SIMD float.
115 * You should normally call the real-precision routine \ref gmx_simd_xor_sign_r.
117 * \param a Values to set sign for
118 * \param b Values used to set sign
119 * \return if b is negative, the sign of a will be changed.
121 * This is equivalent to doing an xor operation on a with the sign bit of b,
122 * with the exception that negative zero is not considered to be negative
123 * on architectures where \ref GMX_SIMD_HAVE_LOGICAL is not set.
125 static gmx_inline gmx_simd_float_t gmx_simdcall
126 gmx_simd_xor_sign_f(gmx_simd_float_t a, gmx_simd_float_t b)
128 #ifdef GMX_SIMD_HAVE_LOGICAL
129 return gmx_simd_xor_f(a, gmx_simd_and_f(gmx_simd_set1_f(GMX_FLOAT_NEGZERO), b));
131 return gmx_simd_blendv_f(a, gmx_simd_fneg_f(a), gmx_simd_cmplt_f(b, gmx_simd_setzero_f()));
135 #ifndef gmx_simd_rsqrt_iter_f
136 /*! \brief Perform one Newton-Raphson iteration to improve 1/sqrt(x) for SIMD float.
138 * This is a low-level routine that should only be used by SIMD math routine
139 * that evaluates the inverse square root.
141 * \param lu Approximation of 1/sqrt(x), typically obtained from lookup.
142 * \param x The reference (starting) value x for which we want 1/sqrt(x).
143 * \return An improved approximation with roughly twice as many bits of accuracy.
145 static gmx_inline gmx_simd_float_t gmx_simdcall
146 gmx_simd_rsqrt_iter_f(gmx_simd_float_t lu, gmx_simd_float_t x)
148 # ifdef GMX_SIMD_HAVE_FMA
149 return gmx_simd_fmadd_f(gmx_simd_fnmadd_f(x, gmx_simd_mul_f(lu, lu), gmx_simd_set1_f(1.0f)), gmx_simd_mul_f(lu, gmx_simd_set1_f(0.5f)), lu);
151 return gmx_simd_mul_f(gmx_simd_set1_f(0.5f), gmx_simd_mul_f(gmx_simd_sub_f(gmx_simd_set1_f(3.0f), gmx_simd_mul_f(gmx_simd_mul_f(lu, lu), x)), lu));
156 /*! \brief Calculate 1/sqrt(x) for SIMD float.
158 * You should normally call the real-precision routine \ref gmx_simd_invsqrt_r.
160 * \param x Argument that must be >0. This routine does not check arguments.
161 * \return 1/sqrt(x). Result is undefined if your argument was invalid.
163 static gmx_inline gmx_simd_float_t gmx_simdcall
164 gmx_simd_invsqrt_f(gmx_simd_float_t x)
166 gmx_simd_float_t lu = gmx_simd_rsqrt_f(x);
167 #if (GMX_SIMD_RSQRT_BITS < GMX_SIMD_MATH_TARGET_SINGLE_BITS)
168 lu = gmx_simd_rsqrt_iter_f(lu, x);
170 #if (GMX_SIMD_RSQRT_BITS*2 < GMX_SIMD_MATH_TARGET_SINGLE_BITS)
171 lu = gmx_simd_rsqrt_iter_f(lu, x);
173 #if (GMX_SIMD_RSQRT_BITS*4 < GMX_SIMD_MATH_TARGET_SINGLE_BITS)
174 lu = gmx_simd_rsqrt_iter_f(lu, x);
179 /*! \brief Calculate 1/sqrt(x) for two SIMD floats.
181 * You should normally call the real-precision routine \ref gmx_simd_invsqrt_pair_r.
183 * \param x0 First set of arguments, x0 must be positive - no argument checking.
184 * \param x1 Second set of arguments, x1 must be positive - no argument checking.
185 * \param[out] out0 Result 1/sqrt(x0)
186 * \param[out] out1 Result 1/sqrt(x1)
188 * In particular for double precision we can sometimes calculate square root
189 * pairs slightly faster by using single precision until the very last step.
191 static gmx_inline void gmx_simdcall
192 gmx_simd_invsqrt_pair_f(gmx_simd_float_t x0, gmx_simd_float_t x1,
193 gmx_simd_float_t *out0, gmx_simd_float_t *out1)
195 *out0 = gmx_simd_invsqrt_f(x0);
196 *out1 = gmx_simd_invsqrt_f(x1);
199 #ifndef gmx_simd_rcp_iter_f
200 /*! \brief Perform one Newton-Raphson iteration to improve 1/x for SIMD float.
202 * This is a low-level routine that should only be used by SIMD math routine
203 * that evaluates the reciprocal.
205 * \param lu Approximation of 1/x, typically obtained from lookup.
206 * \param x The reference (starting) value x for which we want 1/x.
207 * \return An improved approximation with roughly twice as many bits of accuracy.
209 static gmx_inline gmx_simd_float_t gmx_simdcall
210 gmx_simd_rcp_iter_f(gmx_simd_float_t lu, gmx_simd_float_t x)
212 return gmx_simd_mul_f(lu, gmx_simd_fnmadd_f(lu, x, gmx_simd_set1_f(2.0f)));
216 /*! \brief Calculate 1/x for SIMD float.
218 * You should normally call the real-precision routine \ref gmx_simd_inv_r.
220 * \param x Argument that must be nonzero. This routine does not check arguments.
221 * \return 1/x. Result is undefined if your argument was invalid.
223 static gmx_inline gmx_simd_float_t gmx_simdcall
224 gmx_simd_inv_f(gmx_simd_float_t x)
226 gmx_simd_float_t lu = gmx_simd_rcp_f(x);
227 #if (GMX_SIMD_RCP_BITS < GMX_SIMD_MATH_TARGET_SINGLE_BITS)
228 lu = gmx_simd_rcp_iter_f(lu, x);
230 #if (GMX_SIMD_RCP_BITS*2 < GMX_SIMD_MATH_TARGET_SINGLE_BITS)
231 lu = gmx_simd_rcp_iter_f(lu, x);
233 #if (GMX_SIMD_RCP_BITS*4 < GMX_SIMD_MATH_TARGET_SINGLE_BITS)
234 lu = gmx_simd_rcp_iter_f(lu, x);
239 /*! \brief Calculate sqrt(x) correctly for SIMD floats, including argument 0.0.
241 * You should normally call the real-precision routine \ref gmx_simd_sqrt_r.
243 * \param x Argument that must be >=0.
244 * \return sqrt(x). If x=0, the result will correctly be set to 0.
245 * The result is undefined if the input value is negative.
247 static gmx_inline gmx_simd_float_t gmx_simdcall
248 gmx_simd_sqrt_f(gmx_simd_float_t x)
250 gmx_simd_fbool_t mask;
251 gmx_simd_float_t res;
253 mask = gmx_simd_cmpeq_f(x, gmx_simd_setzero_f());
254 res = gmx_simd_blendnotzero_f(gmx_simd_invsqrt_f(x), mask);
255 return gmx_simd_mul_f(res, x);
258 /*! \brief SIMD float log(x). This is the natural logarithm.
260 * You should normally call the real-precision routine \ref gmx_simd_log_r.
262 * \param x Argument, should be >0.
263 * \result The natural logarithm of x. Undefined if argument is invalid.
265 #ifndef gmx_simd_log_f
266 static gmx_inline gmx_simd_float_t gmx_simdcall
267 gmx_simd_log_f(gmx_simd_float_t x)
269 const gmx_simd_float_t half = gmx_simd_set1_f(0.5f);
270 const gmx_simd_float_t one = gmx_simd_set1_f(1.0f);
271 const gmx_simd_float_t sqrt2 = gmx_simd_set1_f(sqrt(2.0f));
272 const gmx_simd_float_t corr = gmx_simd_set1_f(0.693147180559945286226764f);
273 const gmx_simd_float_t CL9 = gmx_simd_set1_f(0.2371599674224853515625f);
274 const gmx_simd_float_t CL7 = gmx_simd_set1_f(0.285279005765914916992188f);
275 const gmx_simd_float_t CL5 = gmx_simd_set1_f(0.400005519390106201171875f);
276 const gmx_simd_float_t CL3 = gmx_simd_set1_f(0.666666567325592041015625f);
277 const gmx_simd_float_t CL1 = gmx_simd_set1_f(2.0f);
278 gmx_simd_float_t fexp, x2, p;
279 gmx_simd_fbool_t mask;
281 fexp = gmx_simd_get_exponent_f(x);
282 x = gmx_simd_get_mantissa_f(x);
284 mask = gmx_simd_cmplt_f(sqrt2, x);
285 /* Adjust to non-IEEE format for x>sqrt(2): exponent += 1, mantissa *= 0.5 */
286 fexp = gmx_simd_add_f(fexp, gmx_simd_blendzero_f(one, mask));
287 x = gmx_simd_mul_f(x, gmx_simd_blendv_f(one, half, mask));
289 x = gmx_simd_mul_f( gmx_simd_sub_f(x, one), gmx_simd_inv_f( gmx_simd_add_f(x, one) ) );
290 x2 = gmx_simd_mul_f(x, x);
292 p = gmx_simd_fmadd_f(CL9, x2, CL7);
293 p = gmx_simd_fmadd_f(p, x2, CL5);
294 p = gmx_simd_fmadd_f(p, x2, CL3);
295 p = gmx_simd_fmadd_f(p, x2, CL1);
296 p = gmx_simd_fmadd_f(p, x, gmx_simd_mul_f(corr, fexp));
302 #ifndef gmx_simd_exp2_f
303 /*! \brief SIMD float 2^x.
305 * You should normally call the real-precision routine \ref gmx_simd_exp2_r.
308 * \result 2^x. Undefined if input argument caused overflow.
310 static gmx_inline gmx_simd_float_t gmx_simdcall
311 gmx_simd_exp2_f(gmx_simd_float_t x)
313 /* Lower bound: Disallow numbers that would lead to an IEEE fp exponent reaching +-127. */
314 const gmx_simd_float_t arglimit = gmx_simd_set1_f(126.0f);
315 const gmx_simd_float_t CC6 = gmx_simd_set1_f(0.0001534581200287996416911311);
316 const gmx_simd_float_t CC5 = gmx_simd_set1_f(0.001339993121934088894618990);
317 const gmx_simd_float_t CC4 = gmx_simd_set1_f(0.009618488957115180159497841);
318 const gmx_simd_float_t CC3 = gmx_simd_set1_f(0.05550328776964726865751735);
319 const gmx_simd_float_t CC2 = gmx_simd_set1_f(0.2402264689063408646490722);
320 const gmx_simd_float_t CC1 = gmx_simd_set1_f(0.6931472057372680777553816);
321 const gmx_simd_float_t one = gmx_simd_set1_f(1.0f);
323 gmx_simd_float_t fexppart;
324 gmx_simd_float_t intpart;
326 gmx_simd_fbool_t valuemask;
328 fexppart = gmx_simd_set_exponent_f(x);
329 intpart = gmx_simd_round_f(x);
330 valuemask = gmx_simd_cmple_f(gmx_simd_fabs_f(x), arglimit);
331 fexppart = gmx_simd_blendzero_f(fexppart, valuemask);
332 x = gmx_simd_sub_f(x, intpart);
334 p = gmx_simd_fmadd_f(CC6, x, CC5);
335 p = gmx_simd_fmadd_f(p, x, CC4);
336 p = gmx_simd_fmadd_f(p, x, CC3);
337 p = gmx_simd_fmadd_f(p, x, CC2);
338 p = gmx_simd_fmadd_f(p, x, CC1);
339 p = gmx_simd_fmadd_f(p, x, one);
340 x = gmx_simd_mul_f(p, fexppart);
345 #ifndef gmx_simd_exp_f
346 /*! \brief SIMD float exp(x).
348 * You should normally call the real-precision routine \ref gmx_simd_exp_r.
350 * In addition to scaling the argument for 2^x this routine correctly does
351 * extended precision arithmetics to improve accuracy.
354 * \result exp(x). Undefined if input argument caused overflow,
355 * which can happen if abs(x) \> 7e13.
357 static gmx_inline gmx_simd_float_t gmx_simdcall
358 gmx_simd_exp_f(gmx_simd_float_t x)
360 const gmx_simd_float_t argscale = gmx_simd_set1_f(1.44269504088896341f);
361 /* Lower bound: Disallow numbers that would lead to an IEEE fp exponent reaching +-127. */
362 const gmx_simd_float_t arglimit = gmx_simd_set1_f(126.0f);
363 const gmx_simd_float_t invargscale0 = gmx_simd_set1_f(0.693145751953125f);
364 const gmx_simd_float_t invargscale1 = gmx_simd_set1_f(1.428606765330187045e-06f);
365 const gmx_simd_float_t CC4 = gmx_simd_set1_f(0.00136324646882712841033936f);
366 const gmx_simd_float_t CC3 = gmx_simd_set1_f(0.00836596917361021041870117f);
367 const gmx_simd_float_t CC2 = gmx_simd_set1_f(0.0416710823774337768554688f);
368 const gmx_simd_float_t CC1 = gmx_simd_set1_f(0.166665524244308471679688f);
369 const gmx_simd_float_t CC0 = gmx_simd_set1_f(0.499999850988388061523438f);
370 const gmx_simd_float_t one = gmx_simd_set1_f(1.0f);
371 gmx_simd_float_t fexppart;
372 gmx_simd_float_t intpart;
373 gmx_simd_float_t y, p;
374 gmx_simd_fbool_t valuemask;
376 y = gmx_simd_mul_f(x, argscale);
377 fexppart = gmx_simd_set_exponent_f(y); /* rounds to nearest int internally */
378 intpart = gmx_simd_round_f(y); /* use same rounding algorithm here */
379 valuemask = gmx_simd_cmple_f(gmx_simd_fabs_f(y), arglimit);
380 fexppart = gmx_simd_blendzero_f(fexppart, valuemask);
382 /* Extended precision arithmetics */
383 x = gmx_simd_fnmadd_f(invargscale0, intpart, x);
384 x = gmx_simd_fnmadd_f(invargscale1, intpart, x);
386 p = gmx_simd_fmadd_f(CC4, x, CC3);
387 p = gmx_simd_fmadd_f(p, x, CC2);
388 p = gmx_simd_fmadd_f(p, x, CC1);
389 p = gmx_simd_fmadd_f(p, x, CC0);
390 p = gmx_simd_fmadd_f(gmx_simd_mul_f(x, x), p, x);
391 p = gmx_simd_add_f(p, one);
392 x = gmx_simd_mul_f(p, fexppart);
397 /*! \brief SIMD float erf(x).
399 * You should normally call the real-precision routine \ref gmx_simd_erf_r.
401 * \param x The value to calculate erf(x) for.
404 * This routine achieves very close to full precision, but we do not care about
405 * the last bit or the subnormal result range.
407 static gmx_inline gmx_simd_float_t gmx_simdcall
408 gmx_simd_erf_f(gmx_simd_float_t x)
410 /* Coefficients for minimax approximation of erf(x)=x*P(x^2) in range [-1,1] */
411 const gmx_simd_float_t CA6 = gmx_simd_set1_f(7.853861353153693e-5f);
412 const gmx_simd_float_t CA5 = gmx_simd_set1_f(-8.010193625184903e-4f);
413 const gmx_simd_float_t CA4 = gmx_simd_set1_f(5.188327685732524e-3f);
414 const gmx_simd_float_t CA3 = gmx_simd_set1_f(-2.685381193529856e-2f);
415 const gmx_simd_float_t CA2 = gmx_simd_set1_f(1.128358514861418e-1f);
416 const gmx_simd_float_t CA1 = gmx_simd_set1_f(-3.761262582423300e-1f);
417 const gmx_simd_float_t CA0 = gmx_simd_set1_f(1.128379165726710f);
418 /* Coefficients for minimax approximation of erfc(x)=Exp(-x^2)*P((1/(x-1))^2) in range [0.67,2] */
419 const gmx_simd_float_t CB9 = gmx_simd_set1_f(-0.0018629930017603923f);
420 const gmx_simd_float_t CB8 = gmx_simd_set1_f(0.003909821287598495f);
421 const gmx_simd_float_t CB7 = gmx_simd_set1_f(-0.0052094582210355615f);
422 const gmx_simd_float_t CB6 = gmx_simd_set1_f(0.005685614362160572f);
423 const gmx_simd_float_t CB5 = gmx_simd_set1_f(-0.0025367682853477272f);
424 const gmx_simd_float_t CB4 = gmx_simd_set1_f(-0.010199799682318782f);
425 const gmx_simd_float_t CB3 = gmx_simd_set1_f(0.04369575504816542f);
426 const gmx_simd_float_t CB2 = gmx_simd_set1_f(-0.11884063474674492f);
427 const gmx_simd_float_t CB1 = gmx_simd_set1_f(0.2732120154030589f);
428 const gmx_simd_float_t CB0 = gmx_simd_set1_f(0.42758357702025784f);
429 /* Coefficients for minimax approximation of erfc(x)=Exp(-x^2)*(1/x)*P((1/x)^2) in range [2,9.19] */
430 const gmx_simd_float_t CC10 = gmx_simd_set1_f(-0.0445555913112064f);
431 const gmx_simd_float_t CC9 = gmx_simd_set1_f(0.21376355144663348f);
432 const gmx_simd_float_t CC8 = gmx_simd_set1_f(-0.3473187200259257f);
433 const gmx_simd_float_t CC7 = gmx_simd_set1_f(0.016690861551248114f);
434 const gmx_simd_float_t CC6 = gmx_simd_set1_f(0.7560973182491192f);
435 const gmx_simd_float_t CC5 = gmx_simd_set1_f(-1.2137903600145787f);
436 const gmx_simd_float_t CC4 = gmx_simd_set1_f(0.8411872321232948f);
437 const gmx_simd_float_t CC3 = gmx_simd_set1_f(-0.08670413896296343f);
438 const gmx_simd_float_t CC2 = gmx_simd_set1_f(-0.27124782687240334f);
439 const gmx_simd_float_t CC1 = gmx_simd_set1_f(-0.0007502488047806069f);
440 const gmx_simd_float_t CC0 = gmx_simd_set1_f(0.5642114853803148f);
441 const gmx_simd_float_t one = gmx_simd_set1_f(1.0f);
442 const gmx_simd_float_t two = gmx_simd_set1_f(2.0f);
444 gmx_simd_float_t x2, x4, y;
445 gmx_simd_float_t t, t2, w, w2;
446 gmx_simd_float_t pA0, pA1, pB0, pB1, pC0, pC1;
447 gmx_simd_float_t expmx2;
448 gmx_simd_float_t res_erf, res_erfc, res;
449 gmx_simd_fbool_t mask;
451 /* Calculate erf() */
452 x2 = gmx_simd_mul_f(x, x);
453 x4 = gmx_simd_mul_f(x2, x2);
455 pA0 = gmx_simd_fmadd_f(CA6, x4, CA4);
456 pA1 = gmx_simd_fmadd_f(CA5, x4, CA3);
457 pA0 = gmx_simd_fmadd_f(pA0, x4, CA2);
458 pA1 = gmx_simd_fmadd_f(pA1, x4, CA1);
459 pA0 = gmx_simd_mul_f(pA0, x4);
460 pA0 = gmx_simd_fmadd_f(pA1, x2, pA0);
461 /* Constant term must come last for precision reasons */
462 pA0 = gmx_simd_add_f(pA0, CA0);
464 res_erf = gmx_simd_mul_f(x, pA0);
467 y = gmx_simd_fabs_f(x);
468 t = gmx_simd_inv_f(y);
469 w = gmx_simd_sub_f(t, one);
470 t2 = gmx_simd_mul_f(t, t);
471 w2 = gmx_simd_mul_f(w, w);
473 /* No need for a floating-point sieve here (as in erfc), since erf()
474 * will never return values that are extremely small for large args.
476 expmx2 = gmx_simd_exp_f( gmx_simd_fneg_f( gmx_simd_mul_f(y, y)));
478 pB1 = gmx_simd_fmadd_f(CB9, w2, CB7);
479 pB0 = gmx_simd_fmadd_f(CB8, w2, CB6);
480 pB1 = gmx_simd_fmadd_f(pB1, w2, CB5);
481 pB0 = gmx_simd_fmadd_f(pB0, w2, CB4);
482 pB1 = gmx_simd_fmadd_f(pB1, w2, CB3);
483 pB0 = gmx_simd_fmadd_f(pB0, w2, CB2);
484 pB1 = gmx_simd_fmadd_f(pB1, w2, CB1);
485 pB0 = gmx_simd_fmadd_f(pB0, w2, CB0);
486 pB0 = gmx_simd_fmadd_f(pB1, w, pB0);
488 pC0 = gmx_simd_fmadd_f(CC10, t2, CC8);
489 pC1 = gmx_simd_fmadd_f(CC9, t2, CC7);
490 pC0 = gmx_simd_fmadd_f(pC0, t2, CC6);
491 pC1 = gmx_simd_fmadd_f(pC1, t2, CC5);
492 pC0 = gmx_simd_fmadd_f(pC0, t2, CC4);
493 pC1 = gmx_simd_fmadd_f(pC1, t2, CC3);
494 pC0 = gmx_simd_fmadd_f(pC0, t2, CC2);
495 pC1 = gmx_simd_fmadd_f(pC1, t2, CC1);
497 pC0 = gmx_simd_fmadd_f(pC0, t2, CC0);
498 pC0 = gmx_simd_fmadd_f(pC1, t, pC0);
499 pC0 = gmx_simd_mul_f(pC0, t);
501 /* SELECT pB0 or pC0 for erfc() */
502 mask = gmx_simd_cmplt_f(two, y);
503 res_erfc = gmx_simd_blendv_f(pB0, pC0, mask);
504 res_erfc = gmx_simd_mul_f(res_erfc, expmx2);
506 /* erfc(x<0) = 2-erfc(|x|) */
507 mask = gmx_simd_cmplt_f(x, gmx_simd_setzero_f());
508 res_erfc = gmx_simd_blendv_f(res_erfc, gmx_simd_sub_f(two, res_erfc), mask);
510 /* Select erf() or erfc() */
511 mask = gmx_simd_cmplt_f(y, gmx_simd_set1_f(0.75f));
512 res = gmx_simd_blendv_f(gmx_simd_sub_f(one, res_erfc), res_erf, mask);
517 /*! \brief SIMD float erfc(x).
519 * You should normally call the real-precision routine \ref gmx_simd_erfc_r.
521 * \param x The value to calculate erfc(x) for.
524 * This routine achieves full precision (bar the last bit) over most of the
525 * input range, but for large arguments where the result is getting close
526 * to the minimum representable numbers we accept slightly larger errors
527 * (think results that are in the ballpark of 10^-30 for single precision,
528 * or 10^-200 for double) since that is not relevant for MD.
530 static gmx_inline gmx_simd_float_t gmx_simdcall
531 gmx_simd_erfc_f(gmx_simd_float_t x)
533 /* Coefficients for minimax approximation of erf(x)=x*P(x^2) in range [-1,1] */
534 const gmx_simd_float_t CA6 = gmx_simd_set1_f(7.853861353153693e-5f);
535 const gmx_simd_float_t CA5 = gmx_simd_set1_f(-8.010193625184903e-4f);
536 const gmx_simd_float_t CA4 = gmx_simd_set1_f(5.188327685732524e-3f);
537 const gmx_simd_float_t CA3 = gmx_simd_set1_f(-2.685381193529856e-2f);
538 const gmx_simd_float_t CA2 = gmx_simd_set1_f(1.128358514861418e-1f);
539 const gmx_simd_float_t CA1 = gmx_simd_set1_f(-3.761262582423300e-1f);
540 const gmx_simd_float_t CA0 = gmx_simd_set1_f(1.128379165726710f);
541 /* Coefficients for minimax approximation of erfc(x)=Exp(-x^2)*P((1/(x-1))^2) in range [0.67,2] */
542 const gmx_simd_float_t CB9 = gmx_simd_set1_f(-0.0018629930017603923f);
543 const gmx_simd_float_t CB8 = gmx_simd_set1_f(0.003909821287598495f);
544 const gmx_simd_float_t CB7 = gmx_simd_set1_f(-0.0052094582210355615f);
545 const gmx_simd_float_t CB6 = gmx_simd_set1_f(0.005685614362160572f);
546 const gmx_simd_float_t CB5 = gmx_simd_set1_f(-0.0025367682853477272f);
547 const gmx_simd_float_t CB4 = gmx_simd_set1_f(-0.010199799682318782f);
548 const gmx_simd_float_t CB3 = gmx_simd_set1_f(0.04369575504816542f);
549 const gmx_simd_float_t CB2 = gmx_simd_set1_f(-0.11884063474674492f);
550 const gmx_simd_float_t CB1 = gmx_simd_set1_f(0.2732120154030589f);
551 const gmx_simd_float_t CB0 = gmx_simd_set1_f(0.42758357702025784f);
552 /* Coefficients for minimax approximation of erfc(x)=Exp(-x^2)*(1/x)*P((1/x)^2) in range [2,9.19] */
553 const gmx_simd_float_t CC10 = gmx_simd_set1_f(-0.0445555913112064f);
554 const gmx_simd_float_t CC9 = gmx_simd_set1_f(0.21376355144663348f);
555 const gmx_simd_float_t CC8 = gmx_simd_set1_f(-0.3473187200259257f);
556 const gmx_simd_float_t CC7 = gmx_simd_set1_f(0.016690861551248114f);
557 const gmx_simd_float_t CC6 = gmx_simd_set1_f(0.7560973182491192f);
558 const gmx_simd_float_t CC5 = gmx_simd_set1_f(-1.2137903600145787f);
559 const gmx_simd_float_t CC4 = gmx_simd_set1_f(0.8411872321232948f);
560 const gmx_simd_float_t CC3 = gmx_simd_set1_f(-0.08670413896296343f);
561 const gmx_simd_float_t CC2 = gmx_simd_set1_f(-0.27124782687240334f);
562 const gmx_simd_float_t CC1 = gmx_simd_set1_f(-0.0007502488047806069f);
563 const gmx_simd_float_t CC0 = gmx_simd_set1_f(0.5642114853803148f);
564 /* Coefficients for expansion of exp(x) in [0,0.1] */
565 /* CD0 and CD1 are both 1.0, so no need to declare them separately */
566 const gmx_simd_float_t CD2 = gmx_simd_set1_f(0.5000066608081202f);
567 const gmx_simd_float_t CD3 = gmx_simd_set1_f(0.1664795422874624f);
568 const gmx_simd_float_t CD4 = gmx_simd_set1_f(0.04379839977652482f);
569 const gmx_simd_float_t one = gmx_simd_set1_f(1.0f);
570 const gmx_simd_float_t two = gmx_simd_set1_f(2.0f);
572 /* We need to use a small trick here, since we cannot assume all SIMD
573 * architectures support integers, and the flag we want (0xfffff000) would
574 * evaluate to NaN (i.e., it cannot be expressed as a floating-point num).
575 * Instead, we represent the flags 0xf0f0f000 and 0x0f0f0000 as valid
576 * fp numbers, and perform a logical or. Since the expression is constant,
577 * we can at least hope it is evaluated at compile-time.
579 #ifdef GMX_SIMD_HAVE_LOGICAL
580 const gmx_simd_float_t sieve = gmx_simd_or_f(gmx_simd_set1_f(-5.965323564e+29f), gmx_simd_set1_f(7.05044434e-30f));
582 const int isieve = 0xFFFFF000;
583 float mem[GMX_SIMD_REAL_WIDTH*2];
584 float * pmem = gmx_simd_align_f(mem);
591 gmx_simd_float_t x2, x4, y;
592 gmx_simd_float_t q, z, t, t2, w, w2;
593 gmx_simd_float_t pA0, pA1, pB0, pB1, pC0, pC1;
594 gmx_simd_float_t expmx2, corr;
595 gmx_simd_float_t res_erf, res_erfc, res;
596 gmx_simd_fbool_t mask;
598 /* Calculate erf() */
599 x2 = gmx_simd_mul_f(x, x);
600 x4 = gmx_simd_mul_f(x2, x2);
602 pA0 = gmx_simd_fmadd_f(CA6, x4, CA4);
603 pA1 = gmx_simd_fmadd_f(CA5, x4, CA3);
604 pA0 = gmx_simd_fmadd_f(pA0, x4, CA2);
605 pA1 = gmx_simd_fmadd_f(pA1, x4, CA1);
606 pA1 = gmx_simd_mul_f(pA1, x2);
607 pA0 = gmx_simd_fmadd_f(pA0, x4, pA1);
608 /* Constant term must come last for precision reasons */
609 pA0 = gmx_simd_add_f(pA0, CA0);
611 res_erf = gmx_simd_mul_f(x, pA0);
614 y = gmx_simd_fabs_f(x);
615 t = gmx_simd_inv_f(y);
616 w = gmx_simd_sub_f(t, one);
617 t2 = gmx_simd_mul_f(t, t);
618 w2 = gmx_simd_mul_f(w, w);
620 * We cannot simply calculate exp(-y2) directly in single precision, since
621 * that will lose a couple of bits of precision due to the multiplication.
622 * Instead, we introduce y=z+w, where the last 12 bits of precision are in w.
623 * Then we get exp(-y2) = exp(-z2)*exp((z-y)*(z+y)).
625 * The only drawback with this is that it requires TWO separate exponential
626 * evaluations, which would be horrible performance-wise. However, the argument
627 * for the second exp() call is always small, so there we simply use a
628 * low-order minimax expansion on [0,0.1].
630 * However, this neat idea requires support for logical ops (and) on
631 * FP numbers, which some vendors decided isn't necessary in their SIMD
632 * instruction sets (Hi, IBM VSX!). In principle we could use some tricks
633 * in double, but we still need memory as a backup when that is not available,
634 * and this case is rare enough that we go directly there...
636 #ifdef GMX_SIMD_HAVE_LOGICAL
637 z = gmx_simd_and_f(y, sieve);
639 gmx_simd_store_f(pmem, y);
640 for (i = 0; i < GMX_SIMD_FLOAT_WIDTH; i++)
643 conv.i = conv.i & isieve;
646 z = gmx_simd_load_f(pmem);
648 q = gmx_simd_mul_f( gmx_simd_sub_f(z, y), gmx_simd_add_f(z, y) );
649 corr = gmx_simd_fmadd_f(CD4, q, CD3);
650 corr = gmx_simd_fmadd_f(corr, q, CD2);
651 corr = gmx_simd_fmadd_f(corr, q, one);
652 corr = gmx_simd_fmadd_f(corr, q, one);
654 expmx2 = gmx_simd_exp_f( gmx_simd_fneg_f( gmx_simd_mul_f(z, z) ) );
655 expmx2 = gmx_simd_mul_f(expmx2, corr);
657 pB1 = gmx_simd_fmadd_f(CB9, w2, CB7);
658 pB0 = gmx_simd_fmadd_f(CB8, w2, CB6);
659 pB1 = gmx_simd_fmadd_f(pB1, w2, CB5);
660 pB0 = gmx_simd_fmadd_f(pB0, w2, CB4);
661 pB1 = gmx_simd_fmadd_f(pB1, w2, CB3);
662 pB0 = gmx_simd_fmadd_f(pB0, w2, CB2);
663 pB1 = gmx_simd_fmadd_f(pB1, w2, CB1);
664 pB0 = gmx_simd_fmadd_f(pB0, w2, CB0);
665 pB0 = gmx_simd_fmadd_f(pB1, w, pB0);
667 pC0 = gmx_simd_fmadd_f(CC10, t2, CC8);
668 pC1 = gmx_simd_fmadd_f(CC9, t2, CC7);
669 pC0 = gmx_simd_fmadd_f(pC0, t2, CC6);
670 pC1 = gmx_simd_fmadd_f(pC1, t2, CC5);
671 pC0 = gmx_simd_fmadd_f(pC0, t2, CC4);
672 pC1 = gmx_simd_fmadd_f(pC1, t2, CC3);
673 pC0 = gmx_simd_fmadd_f(pC0, t2, CC2);
674 pC1 = gmx_simd_fmadd_f(pC1, t2, CC1);
676 pC0 = gmx_simd_fmadd_f(pC0, t2, CC0);
677 pC0 = gmx_simd_fmadd_f(pC1, t, pC0);
678 pC0 = gmx_simd_mul_f(pC0, t);
680 /* SELECT pB0 or pC0 for erfc() */
681 mask = gmx_simd_cmplt_f(two, y);
682 res_erfc = gmx_simd_blendv_f(pB0, pC0, mask);
683 res_erfc = gmx_simd_mul_f(res_erfc, expmx2);
685 /* erfc(x<0) = 2-erfc(|x|) */
686 mask = gmx_simd_cmplt_f(x, gmx_simd_setzero_f());
687 res_erfc = gmx_simd_blendv_f(res_erfc, gmx_simd_sub_f(two, res_erfc), mask);
689 /* Select erf() or erfc() */
690 mask = gmx_simd_cmplt_f(y, gmx_simd_set1_f(0.75f));
691 res = gmx_simd_blendv_f(res_erfc, gmx_simd_sub_f(one, res_erf), mask);
696 /*! \brief SIMD float sin \& cos.
698 * You should normally call the real-precision routine \ref gmx_simd_sincos_r.
700 * \param x The argument to evaluate sin/cos for
701 * \param[out] sinval Sin(x)
702 * \param[out] cosval Cos(x)
704 * This version achieves close to machine precision, but for very large
705 * magnitudes of the argument we inherently begin to lose accuracy due to the
706 * argument reduction, despite using extended precision arithmetics internally.
708 static gmx_inline void gmx_simdcall
709 gmx_simd_sincos_f(gmx_simd_float_t x, gmx_simd_float_t *sinval, gmx_simd_float_t *cosval)
711 /* Constants to subtract Pi/4*x from y while minimizing precision loss */
712 const gmx_simd_float_t argred0 = gmx_simd_set1_f(1.5703125);
713 const gmx_simd_float_t argred1 = gmx_simd_set1_f(4.83751296997070312500e-04f);
714 const gmx_simd_float_t argred2 = gmx_simd_set1_f(7.54953362047672271729e-08f);
715 const gmx_simd_float_t argred3 = gmx_simd_set1_f(2.56334406825708960298e-12f);
716 const gmx_simd_float_t two_over_pi = gmx_simd_set1_f(2.0f/M_PI);
717 const gmx_simd_float_t const_sin2 = gmx_simd_set1_f(-1.9515295891e-4f);
718 const gmx_simd_float_t const_sin1 = gmx_simd_set1_f( 8.3321608736e-3f);
719 const gmx_simd_float_t const_sin0 = gmx_simd_set1_f(-1.6666654611e-1f);
720 const gmx_simd_float_t const_cos2 = gmx_simd_set1_f( 2.443315711809948e-5f);
721 const gmx_simd_float_t const_cos1 = gmx_simd_set1_f(-1.388731625493765e-3f);
722 const gmx_simd_float_t const_cos0 = gmx_simd_set1_f( 4.166664568298827e-2f);
723 const gmx_simd_float_t half = gmx_simd_set1_f(0.5f);
724 const gmx_simd_float_t one = gmx_simd_set1_f(1.0f);
725 gmx_simd_float_t ssign, csign;
726 gmx_simd_float_t x2, y, z, psin, pcos, sss, ccc;
727 gmx_simd_fbool_t mask;
728 #if (defined GMX_SIMD_HAVE_FINT32) && (defined GMX_SIMD_HAVE_FINT32_ARITHMETICS) && (defined GMX_SIMD_HAVE_LOGICAL)
729 const gmx_simd_fint32_t ione = gmx_simd_set1_fi(1);
730 const gmx_simd_fint32_t itwo = gmx_simd_set1_fi(2);
731 gmx_simd_fint32_t iy;
733 z = gmx_simd_mul_f(x, two_over_pi);
734 iy = gmx_simd_cvt_f2i(z);
735 y = gmx_simd_round_f(z);
737 mask = gmx_simd_cvt_fib2fb(gmx_simd_cmpeq_fi(gmx_simd_and_fi(iy, ione), gmx_simd_setzero_fi()));
738 ssign = gmx_simd_blendzero_f(gmx_simd_set1_f(GMX_FLOAT_NEGZERO), gmx_simd_cvt_fib2fb(gmx_simd_cmpeq_fi(gmx_simd_and_fi(iy, itwo), itwo)));
739 csign = gmx_simd_blendzero_f(gmx_simd_set1_f(GMX_FLOAT_NEGZERO), gmx_simd_cvt_fib2fb(gmx_simd_cmpeq_fi(gmx_simd_and_fi(gmx_simd_add_fi(iy, ione), itwo), itwo)));
741 const gmx_simd_float_t quarter = gmx_simd_set1_f(0.25f);
742 const gmx_simd_float_t minusquarter = gmx_simd_set1_f(-0.25f);
744 gmx_simd_fbool_t m1, m2, m3;
746 /* The most obvious way to find the arguments quadrant in the unit circle
747 * to calculate the sign is to use integer arithmetic, but that is not
748 * present in all SIMD implementations. As an alternative, we have devised a
749 * pure floating-point algorithm that uses truncation for argument reduction
750 * so that we get a new value 0<=q<1 over the unit circle, and then
751 * do floating-point comparisons with fractions. This is likely to be
752 * slightly slower (~10%) due to the longer latencies of floating-point, so
753 * we only use it when integer SIMD arithmetic is not present.
756 x = gmx_simd_fabs_f(x);
757 /* It is critical that half-way cases are rounded down */
758 z = gmx_simd_fmadd_f(x, two_over_pi, half);
759 y = gmx_simd_trunc_f(z);
760 q = gmx_simd_mul_f(z, quarter);
761 q = gmx_simd_sub_f(q, gmx_simd_trunc_f(q));
762 /* z now starts at 0.0 for x=-pi/4 (although neg. values cannot occur), and
763 * then increased by 1.0 as x increases by 2*Pi, when it resets to 0.0.
764 * This removes the 2*Pi periodicity without using any integer arithmetic.
765 * First check if y had the value 2 or 3, set csign if true.
767 q = gmx_simd_sub_f(q, half);
768 /* If we have logical operations we can work directly on the signbit, which
769 * saves instructions. Otherwise we need to represent signs as +1.0/-1.0.
770 * Thus, if you are altering defines to debug alternative code paths, the
771 * two GMX_SIMD_HAVE_LOGICAL sections in this routine must either both be
772 * active or inactive - you will get errors if only one is used.
774 # ifdef GMX_SIMD_HAVE_LOGICAL
775 ssign = gmx_simd_and_f(ssign, gmx_simd_set1_f(GMX_FLOAT_NEGZERO));
776 csign = gmx_simd_andnot_f(q, gmx_simd_set1_f(GMX_FLOAT_NEGZERO));
777 ssign = gmx_simd_xor_f(ssign, csign);
779 csign = gmx_simd_xor_sign_f(gmx_simd_set1_f(-1.0f), q);
780 // ALT: csign = gmx_simd_fneg_f(gmx_simd_copysign(gmx_simd_set1_f(1.0),q));
782 ssign = gmx_simd_xor_sign_f(ssign, csign); /* swap ssign if csign was set. */
784 /* Check if y had value 1 or 3 (remember we subtracted 0.5 from q) */
785 m1 = gmx_simd_cmplt_f(q, minusquarter);
786 m2 = gmx_simd_cmple_f(gmx_simd_setzero_f(), q);
787 m3 = gmx_simd_cmplt_f(q, quarter);
788 m2 = gmx_simd_and_fb(m2, m3);
789 mask = gmx_simd_or_fb(m1, m2);
790 /* where mask is FALSE, set sign. */
791 csign = gmx_simd_xor_sign_f(csign, gmx_simd_blendv_f(gmx_simd_set1_f(-1.0f), one, mask));
793 x = gmx_simd_fnmadd_f(y, argred0, x);
794 x = gmx_simd_fnmadd_f(y, argred1, x);
795 x = gmx_simd_fnmadd_f(y, argred2, x);
796 x = gmx_simd_fnmadd_f(y, argred3, x);
797 x2 = gmx_simd_mul_f(x, x);
799 psin = gmx_simd_fmadd_f(const_sin2, x2, const_sin1);
800 psin = gmx_simd_fmadd_f(psin, x2, const_sin0);
801 psin = gmx_simd_fmadd_f(psin, gmx_simd_mul_f(x, x2), x);
802 pcos = gmx_simd_fmadd_f(const_cos2, x2, const_cos1);
803 pcos = gmx_simd_fmadd_f(pcos, x2, const_cos0);
804 pcos = gmx_simd_fmsub_f(pcos, x2, half);
805 pcos = gmx_simd_fmadd_f(pcos, x2, one);
807 sss = gmx_simd_blendv_f(pcos, psin, mask);
808 ccc = gmx_simd_blendv_f(psin, pcos, mask);
809 /* See comment for GMX_SIMD_HAVE_LOGICAL section above. */
810 #ifdef GMX_SIMD_HAVE_LOGICAL
811 *sinval = gmx_simd_xor_f(sss, ssign);
812 *cosval = gmx_simd_xor_f(ccc, csign);
814 *sinval = gmx_simd_xor_sign_f(sss, ssign);
815 *cosval = gmx_simd_xor_sign_f(ccc, csign);
819 /*! \brief SIMD float sin(x).
821 * You should normally call the real-precision routine \ref gmx_simd_sin_r.
823 * \param x The argument to evaluate sin for
826 * \attention Do NOT call both sin & cos if you need both results, since each of them
827 * will then call \ref gmx_simd_sincos_r and waste a factor 2 in performance.
829 static gmx_inline gmx_simd_float_t gmx_simdcall
830 gmx_simd_sin_f(gmx_simd_float_t x)
832 gmx_simd_float_t s, c;
833 gmx_simd_sincos_f(x, &s, &c);
837 /*! \brief SIMD float cos(x).
839 * You should normally call the real-precision routine \ref gmx_simd_cos_r.
841 * \param x The argument to evaluate cos for
844 * \attention Do NOT call both sin & cos if you need both results, since each of them
845 * will then call \ref gmx_simd_sincos_r and waste a factor 2 in performance.
847 static gmx_inline gmx_simd_float_t gmx_simdcall
848 gmx_simd_cos_f(gmx_simd_float_t x)
850 gmx_simd_float_t s, c;
851 gmx_simd_sincos_f(x, &s, &c);
855 /*! \brief SIMD float tan(x).
857 * You should normally call the real-precision routine \ref gmx_simd_tan_r.
859 * \param x The argument to evaluate tan for
862 static gmx_inline gmx_simd_float_t gmx_simdcall
863 gmx_simd_tan_f(gmx_simd_float_t x)
865 const gmx_simd_float_t argred0 = gmx_simd_set1_f(1.5703125);
866 const gmx_simd_float_t argred1 = gmx_simd_set1_f(4.83751296997070312500e-04f);
867 const gmx_simd_float_t argred2 = gmx_simd_set1_f(7.54953362047672271729e-08f);
868 const gmx_simd_float_t argred3 = gmx_simd_set1_f(2.56334406825708960298e-12f);
869 const gmx_simd_float_t two_over_pi = gmx_simd_set1_f(2.0f/M_PI);
870 const gmx_simd_float_t CT6 = gmx_simd_set1_f(0.009498288995810566122993911);
871 const gmx_simd_float_t CT5 = gmx_simd_set1_f(0.002895755790837379295226923);
872 const gmx_simd_float_t CT4 = gmx_simd_set1_f(0.02460087336161924491836265);
873 const gmx_simd_float_t CT3 = gmx_simd_set1_f(0.05334912882656359828045988);
874 const gmx_simd_float_t CT2 = gmx_simd_set1_f(0.1333989091464957704418495);
875 const gmx_simd_float_t CT1 = gmx_simd_set1_f(0.3333307599244198227797507);
877 gmx_simd_float_t x2, p, y, z;
878 gmx_simd_fbool_t mask;
880 #if (defined GMX_SIMD_HAVE_FINT32) && (defined GMX_SIMD_HAVE_FINT32_ARITHMETICS) && (defined GMX_SIMD_HAVE_LOGICAL)
881 gmx_simd_fint32_t iy;
882 gmx_simd_fint32_t ione = gmx_simd_set1_fi(1);
884 z = gmx_simd_mul_f(x, two_over_pi);
885 iy = gmx_simd_cvt_f2i(z);
886 y = gmx_simd_round_f(z);
887 mask = gmx_simd_cvt_fib2fb(gmx_simd_cmpeq_fi(gmx_simd_and_fi(iy, ione), ione));
889 x = gmx_simd_fnmadd_f(y, argred0, x);
890 x = gmx_simd_fnmadd_f(y, argred1, x);
891 x = gmx_simd_fnmadd_f(y, argred2, x);
892 x = gmx_simd_fnmadd_f(y, argred3, x);
893 x = gmx_simd_xor_f(gmx_simd_blendzero_f(gmx_simd_set1_f(GMX_FLOAT_NEGZERO), mask), x);
895 const gmx_simd_float_t quarter = gmx_simd_set1_f(0.25f);
896 const gmx_simd_float_t half = gmx_simd_set1_f(0.5f);
897 const gmx_simd_float_t threequarter = gmx_simd_set1_f(0.75f);
898 gmx_simd_float_t w, q;
899 gmx_simd_fbool_t m1, m2, m3;
901 w = gmx_simd_fabs_f(x);
902 z = gmx_simd_fmadd_f(w, two_over_pi, half);
903 y = gmx_simd_trunc_f(z);
904 q = gmx_simd_mul_f(z, quarter);
905 q = gmx_simd_sub_f(q, gmx_simd_trunc_f(q));
906 m1 = gmx_simd_cmple_f(quarter, q);
907 m2 = gmx_simd_cmplt_f(q, half);
908 m3 = gmx_simd_cmple_f(threequarter, q);
909 m1 = gmx_simd_and_fb(m1, m2);
910 mask = gmx_simd_or_fb(m1, m3);
911 w = gmx_simd_fnmadd_f(y, argred0, w);
912 w = gmx_simd_fnmadd_f(y, argred1, w);
913 w = gmx_simd_fnmadd_f(y, argred2, w);
914 w = gmx_simd_fnmadd_f(y, argred3, w);
916 w = gmx_simd_blendv_f(w, gmx_simd_fneg_f(w), mask);
917 x = gmx_simd_xor_sign_f(w, x);
919 x2 = gmx_simd_mul_f(x, x);
920 p = gmx_simd_fmadd_f(CT6, x2, CT5);
921 p = gmx_simd_fmadd_f(p, x2, CT4);
922 p = gmx_simd_fmadd_f(p, x2, CT3);
923 p = gmx_simd_fmadd_f(p, x2, CT2);
924 p = gmx_simd_fmadd_f(p, x2, CT1);
925 p = gmx_simd_fmadd_f(x2, gmx_simd_mul_f(p, x), x);
927 p = gmx_simd_blendv_f( p, gmx_simd_inv_f(p), mask);
931 /*! \brief SIMD float asin(x).
933 * You should normally call the real-precision routine \ref gmx_simd_asin_r.
935 * \param x The argument to evaluate asin for
938 static gmx_inline gmx_simd_float_t gmx_simdcall
939 gmx_simd_asin_f(gmx_simd_float_t x)
941 const gmx_simd_float_t limitlow = gmx_simd_set1_f(1e-4f);
942 const gmx_simd_float_t half = gmx_simd_set1_f(0.5f);
943 const gmx_simd_float_t one = gmx_simd_set1_f(1.0f);
944 const gmx_simd_float_t halfpi = gmx_simd_set1_f((float)M_PI/2.0f);
945 const gmx_simd_float_t CC5 = gmx_simd_set1_f(4.2163199048E-2f);
946 const gmx_simd_float_t CC4 = gmx_simd_set1_f(2.4181311049E-2f);
947 const gmx_simd_float_t CC3 = gmx_simd_set1_f(4.5470025998E-2f);
948 const gmx_simd_float_t CC2 = gmx_simd_set1_f(7.4953002686E-2f);
949 const gmx_simd_float_t CC1 = gmx_simd_set1_f(1.6666752422E-1f);
950 gmx_simd_float_t xabs;
951 gmx_simd_float_t z, z1, z2, q, q1, q2;
952 gmx_simd_float_t pA, pB;
953 gmx_simd_fbool_t mask;
955 xabs = gmx_simd_fabs_f(x);
956 mask = gmx_simd_cmplt_f(half, xabs);
957 z1 = gmx_simd_mul_f(half, gmx_simd_sub_f(one, xabs));
958 q1 = gmx_simd_mul_f(z1, gmx_simd_invsqrt_f(z1));
959 q1 = gmx_simd_blendnotzero_f(q1, gmx_simd_cmpeq_f(xabs, one));
961 z2 = gmx_simd_mul_f(q2, q2);
962 z = gmx_simd_blendv_f(z2, z1, mask);
963 q = gmx_simd_blendv_f(q2, q1, mask);
965 z2 = gmx_simd_mul_f(z, z);
966 pA = gmx_simd_fmadd_f(CC5, z2, CC3);
967 pB = gmx_simd_fmadd_f(CC4, z2, CC2);
968 pA = gmx_simd_fmadd_f(pA, z2, CC1);
969 pA = gmx_simd_mul_f(pA, z);
970 z = gmx_simd_fmadd_f(pB, z2, pA);
971 z = gmx_simd_fmadd_f(z, q, q);
972 q2 = gmx_simd_sub_f(halfpi, z);
973 q2 = gmx_simd_sub_f(q2, z);
974 z = gmx_simd_blendv_f(z, q2, mask);
976 mask = gmx_simd_cmplt_f(limitlow, xabs);
977 z = gmx_simd_blendv_f( xabs, z, mask );
978 z = gmx_simd_xor_sign_f(z, x);
983 /*! \brief SIMD float acos(x).
985 * You should normally call the real-precision routine \ref gmx_simd_acos_r.
987 * \param x The argument to evaluate acos for
990 static gmx_inline gmx_simd_float_t gmx_simdcall
991 gmx_simd_acos_f(gmx_simd_float_t x)
993 const gmx_simd_float_t one = gmx_simd_set1_f(1.0f);
994 const gmx_simd_float_t half = gmx_simd_set1_f(0.5f);
995 const gmx_simd_float_t pi = gmx_simd_set1_f((float)M_PI);
996 const gmx_simd_float_t halfpi = gmx_simd_set1_f((float)M_PI/2.0f);
997 gmx_simd_float_t xabs;
998 gmx_simd_float_t z, z1, z2, z3;
999 gmx_simd_fbool_t mask1, mask2;
1001 xabs = gmx_simd_fabs_f(x);
1002 mask1 = gmx_simd_cmplt_f(half, xabs);
1003 mask2 = gmx_simd_cmplt_f(gmx_simd_setzero_f(), x);
1005 z = gmx_simd_mul_f(half, gmx_simd_sub_f(one, xabs));
1006 z = gmx_simd_mul_f(z, gmx_simd_invsqrt_f(z));
1007 z = gmx_simd_blendnotzero_f(z, gmx_simd_cmpeq_f(xabs, one));
1008 z = gmx_simd_blendv_f(x, z, mask1);
1009 z = gmx_simd_asin_f(z);
1011 z2 = gmx_simd_add_f(z, z);
1012 z1 = gmx_simd_sub_f(pi, z2);
1013 z3 = gmx_simd_sub_f(halfpi, z);
1014 z = gmx_simd_blendv_f(z1, z2, mask2);
1015 z = gmx_simd_blendv_f(z3, z, mask1);
1020 /*! \brief SIMD float asin(x).
1022 * You should normally call the real-precision routine \ref gmx_simd_atan_r.
1024 * \param x The argument to evaluate atan for
1025 * \result Atan(x), same argument/value range as standard math library.
1027 static gmx_inline gmx_simd_float_t gmx_simdcall
1028 gmx_simd_atan_f(gmx_simd_float_t x)
1030 const gmx_simd_float_t halfpi = gmx_simd_set1_f(M_PI/2);
1031 const gmx_simd_float_t CA17 = gmx_simd_set1_f(0.002823638962581753730774f);
1032 const gmx_simd_float_t CA15 = gmx_simd_set1_f(-0.01595690287649631500244f);
1033 const gmx_simd_float_t CA13 = gmx_simd_set1_f(0.04250498861074447631836f);
1034 const gmx_simd_float_t CA11 = gmx_simd_set1_f(-0.07489009201526641845703f);
1035 const gmx_simd_float_t CA9 = gmx_simd_set1_f(0.1063479334115982055664f);
1036 const gmx_simd_float_t CA7 = gmx_simd_set1_f(-0.1420273631811141967773f);
1037 const gmx_simd_float_t CA5 = gmx_simd_set1_f(0.1999269574880599975585f);
1038 const gmx_simd_float_t CA3 = gmx_simd_set1_f(-0.3333310186862945556640f);
1039 gmx_simd_float_t x2, x3, x4, pA, pB;
1040 gmx_simd_fbool_t mask, mask2;
1042 mask = gmx_simd_cmplt_f(x, gmx_simd_setzero_f());
1043 x = gmx_simd_fabs_f(x);
1044 mask2 = gmx_simd_cmplt_f(gmx_simd_set1_f(1.0f), x);
1045 x = gmx_simd_blendv_f(x, gmx_simd_inv_f(x), mask2);
1047 x2 = gmx_simd_mul_f(x, x);
1048 x3 = gmx_simd_mul_f(x2, x);
1049 x4 = gmx_simd_mul_f(x2, x2);
1050 pA = gmx_simd_fmadd_f(CA17, x4, CA13);
1051 pB = gmx_simd_fmadd_f(CA15, x4, CA11);
1052 pA = gmx_simd_fmadd_f(pA, x4, CA9);
1053 pB = gmx_simd_fmadd_f(pB, x4, CA7);
1054 pA = gmx_simd_fmadd_f(pA, x4, CA5);
1055 pB = gmx_simd_fmadd_f(pB, x4, CA3);
1056 pA = gmx_simd_fmadd_f(pA, x2, pB);
1057 pA = gmx_simd_fmadd_f(pA, x3, x);
1059 pA = gmx_simd_blendv_f(pA, gmx_simd_sub_f(halfpi, pA), mask2);
1060 pA = gmx_simd_blendv_f(pA, gmx_simd_fneg_f(pA), mask);
1065 /*! \brief SIMD float atan2(y,x).
1067 * You should normally call the real-precision routine \ref gmx_simd_atan2_r.
1069 * \param y Y component of vector, any quartile
1070 * \param x X component of vector, any quartile
1071 * \result Atan(y,x), same argument/value range as standard math library.
1073 * \note This routine should provide correct results for all finite
1074 * non-zero or positive-zero arguments. However, negative zero arguments will
1075 * be treated as positive zero, which means the return value will deviate from
1076 * the standard math library atan2(y,x) for those cases. That should not be
1077 * of any concern in Gromacs, and in particular it will not affect calculations
1078 * of angles from vectors.
1080 static gmx_inline gmx_simd_float_t gmx_simdcall
1081 gmx_simd_atan2_f(gmx_simd_float_t y, gmx_simd_float_t x)
1083 const gmx_simd_float_t pi = gmx_simd_set1_f(M_PI);
1084 const gmx_simd_float_t halfpi = gmx_simd_set1_f(M_PI/2.0);
1085 gmx_simd_float_t xinv, p, aoffset;
1086 gmx_simd_fbool_t mask_x0, mask_y0, mask_xlt0, mask_ylt0;
1088 mask_x0 = gmx_simd_cmpeq_f(x, gmx_simd_setzero_f());
1089 mask_y0 = gmx_simd_cmpeq_f(y, gmx_simd_setzero_f());
1090 mask_xlt0 = gmx_simd_cmplt_f(x, gmx_simd_setzero_f());
1091 mask_ylt0 = gmx_simd_cmplt_f(y, gmx_simd_setzero_f());
1093 aoffset = gmx_simd_blendzero_f(halfpi, mask_x0);
1094 aoffset = gmx_simd_blendnotzero_f(aoffset, mask_y0);
1096 aoffset = gmx_simd_blendv_f(aoffset, pi, mask_xlt0);
1097 aoffset = gmx_simd_blendv_f(aoffset, gmx_simd_fneg_f(aoffset), mask_ylt0);
1099 xinv = gmx_simd_blendnotzero_f(gmx_simd_inv_f(x), mask_x0);
1100 p = gmx_simd_mul_f(y, xinv);
1101 p = gmx_simd_atan_f(p);
1102 p = gmx_simd_add_f(p, aoffset);
1107 /*! \brief Calculate the force correction due to PME analytically in SIMD float.
1109 * You should normally call the real-precision routine \ref gmx_simd_pmecorrF_r.
1111 * \param z2 \f$(r \beta)^2\f$ - see below for details.
1112 * \result Correction factor to coulomb force - see below for details.
1114 * This routine is meant to enable analytical evaluation of the
1115 * direct-space PME electrostatic force to avoid tables.
1117 * The direct-space potential should be \f$ \mbox{erfc}(\beta r)/r\f$, but there
1118 * are some problems evaluating that:
1120 * First, the error function is difficult (read: expensive) to
1121 * approxmiate accurately for intermediate to large arguments, and
1122 * this happens already in ranges of \f$(\beta r)\f$ that occur in simulations.
1123 * Second, we now try to avoid calculating potentials in Gromacs but
1124 * use forces directly.
1126 * We can simply things slight by noting that the PME part is really
1127 * a correction to the normal Coulomb force since \f$\mbox{erfc}(z)=1-\mbox{erf}(z)\f$, i.e.
1129 * V = \frac{1}{r} - \frac{\mbox{erf}(\beta r)}{r}
1131 * The first term we already have from the inverse square root, so
1132 * that we can leave out of this routine.
1134 * For pme tolerances of 1e-3 to 1e-8 and cutoffs of 0.5nm to 1.8nm,
1135 * the argument \f$beta r\f$ will be in the range 0.15 to ~4, which is
1136 * the range used for the minimax fit. Use your favorite plotting program
1137 * to realize how well-behaved \f$\frac{\mbox{erf}(z)}{z}\f$ is in this range!
1139 * We approximate \f$f(z)=\mbox{erf}(z)/z\f$ with a rational minimax polynomial.
1140 * However, it turns out it is more efficient to approximate \f$f(z)/z\f$ and
1141 * then only use even powers. This is another minor optimization, since
1142 * we actually \a want \f$f(z)/z\f$, because it is going to be multiplied by
1143 * the vector between the two atoms to get the vectorial force. The
1144 * fastest flops are the ones we can avoid calculating!
1146 * So, here's how it should be used:
1148 * 1. Calculate \f$r^2\f$.
1149 * 2. Multiply by \f$\beta^2\f$, so you get \f$z^2=(\beta r)^2\f$.
1150 * 3. Evaluate this routine with \f$z^2\f$ as the argument.
1151 * 4. The return value is the expression:
1154 * \frac{2 \exp{-z^2}}{\sqrt{\pi} z^2}-\frac{\mbox{erf}(z)}{z^3}
1157 * 5. Multiply the entire expression by \f$\beta^3\f$. This will get you
1160 * \frac{2 \beta^3 \exp(-z^2)}{\sqrt{\pi} z^2} - \frac{\beta^3 \mbox{erf}(z)}{z^3}
1163 * or, switching back to \f$r\f$ (since \f$z=r \beta\f$):
1166 * \frac{2 \beta \exp(-r^2 \beta^2)}{\sqrt{\pi} r^2} - \frac{\mbox{erf}(r \beta)}{r^3}
1169 * With a bit of math exercise you should be able to confirm that
1173 * \frac{\frac{d}{dr}\left( \frac{\mbox{erf}(\beta r)}{r} \right)}{r}
1176 * 6. Add the result to \f$r^{-3}\f$, multiply by the product of the charges,
1177 * and you have your force (divided by \f$r\f$). A final multiplication
1178 * with the vector connecting the two particles and you have your
1179 * vectorial force to add to the particles.
1181 * This approximation achieves an error slightly lower than 1e-6
1182 * in single precision and 1e-11 in double precision
1183 * for arguments smaller than 16 (\f$\beta r \leq 4 \f$);
1184 * when added to \f$1/r\f$ the error will be insignificant.
1185 * For \f$\beta r \geq 7206\f$ the return value can be inf or NaN.
1188 static gmx_inline gmx_simd_float_t gmx_simdcall
1189 gmx_simd_pmecorrF_f(gmx_simd_float_t z2)
1191 const gmx_simd_float_t FN6 = gmx_simd_set1_f(-1.7357322914161492954e-8f);
1192 const gmx_simd_float_t FN5 = gmx_simd_set1_f(1.4703624142580877519e-6f);
1193 const gmx_simd_float_t FN4 = gmx_simd_set1_f(-0.000053401640219807709149f);
1194 const gmx_simd_float_t FN3 = gmx_simd_set1_f(0.0010054721316683106153f);
1195 const gmx_simd_float_t FN2 = gmx_simd_set1_f(-0.019278317264888380590f);
1196 const gmx_simd_float_t FN1 = gmx_simd_set1_f(0.069670166153766424023f);
1197 const gmx_simd_float_t FN0 = gmx_simd_set1_f(-0.75225204789749321333f);
1199 const gmx_simd_float_t FD4 = gmx_simd_set1_f(0.0011193462567257629232f);
1200 const gmx_simd_float_t FD3 = gmx_simd_set1_f(0.014866955030185295499f);
1201 const gmx_simd_float_t FD2 = gmx_simd_set1_f(0.11583842382862377919f);
1202 const gmx_simd_float_t FD1 = gmx_simd_set1_f(0.50736591960530292870f);
1203 const gmx_simd_float_t FD0 = gmx_simd_set1_f(1.0f);
1205 gmx_simd_float_t z4;
1206 gmx_simd_float_t polyFN0, polyFN1, polyFD0, polyFD1;
1208 z4 = gmx_simd_mul_f(z2, z2);
1210 polyFD0 = gmx_simd_fmadd_f(FD4, z4, FD2);
1211 polyFD1 = gmx_simd_fmadd_f(FD3, z4, FD1);
1212 polyFD0 = gmx_simd_fmadd_f(polyFD0, z4, FD0);
1213 polyFD0 = gmx_simd_fmadd_f(polyFD1, z2, polyFD0);
1215 polyFD0 = gmx_simd_inv_f(polyFD0);
1217 polyFN0 = gmx_simd_fmadd_f(FN6, z4, FN4);
1218 polyFN1 = gmx_simd_fmadd_f(FN5, z4, FN3);
1219 polyFN0 = gmx_simd_fmadd_f(polyFN0, z4, FN2);
1220 polyFN1 = gmx_simd_fmadd_f(polyFN1, z4, FN1);
1221 polyFN0 = gmx_simd_fmadd_f(polyFN0, z4, FN0);
1222 polyFN0 = gmx_simd_fmadd_f(polyFN1, z2, polyFN0);
1224 return gmx_simd_mul_f(polyFN0, polyFD0);
1229 /*! \brief Calculate the potential correction due to PME analytically in SIMD float.
1231 * You should normally call the real-precision routine \ref gmx_simd_pmecorrV_r.
1233 * \param z2 \f$(r \beta)^2\f$ - see below for details.
1234 * \result Correction factor to coulomb potential - see below for details.
1236 * See \ref gmx_simd_pmecorrF_f for details about the approximation.
1238 * This routine calculates \f$\mbox{erf}(z)/z\f$, although you should provide \f$z^2\f$
1239 * as the input argument.
1241 * Here's how it should be used:
1243 * 1. Calculate \f$r^2\f$.
1244 * 2. Multiply by \f$\beta^2\f$, so you get \f$z^2=\beta^2*r^2\f$.
1245 * 3. Evaluate this routine with z^2 as the argument.
1246 * 4. The return value is the expression:
1249 * \frac{\mbox{erf}(z)}{z}
1252 * 5. Multiply the entire expression by beta and switching back to \f$r\f$ (since \f$z=r \beta\f$):
1255 * \frac{\mbox{erf}(r \beta)}{r}
1258 * 6. Subtract the result from \f$1/r\f$, multiply by the product of the charges,
1259 * and you have your potential.
1261 * This approximation achieves an error slightly lower than 1e-6
1262 * in single precision and 4e-11 in double precision
1263 * for arguments smaller than 16 (\f$ 0.15 \leq \beta r \leq 4 \f$);
1264 * for \f$ \beta r \leq 0.15\f$ the error can be twice as high;
1265 * when added to \f$1/r\f$ the error will be insignificant.
1266 * For \f$\beta r \geq 7142\f$ the return value can be inf or NaN.
1268 static gmx_inline gmx_simd_float_t gmx_simdcall
1269 gmx_simd_pmecorrV_f(gmx_simd_float_t z2)
1271 const gmx_simd_float_t VN6 = gmx_simd_set1_f(1.9296833005951166339e-8f);
1272 const gmx_simd_float_t VN5 = gmx_simd_set1_f(-1.4213390571557850962e-6f);
1273 const gmx_simd_float_t VN4 = gmx_simd_set1_f(0.000041603292906656984871f);
1274 const gmx_simd_float_t VN3 = gmx_simd_set1_f(-0.00013134036773265025626f);
1275 const gmx_simd_float_t VN2 = gmx_simd_set1_f(0.038657983986041781264f);
1276 const gmx_simd_float_t VN1 = gmx_simd_set1_f(0.11285044772717598220f);
1277 const gmx_simd_float_t VN0 = gmx_simd_set1_f(1.1283802385263030286f);
1279 const gmx_simd_float_t VD3 = gmx_simd_set1_f(0.0066752224023576045451f);
1280 const gmx_simd_float_t VD2 = gmx_simd_set1_f(0.078647795836373922256f);
1281 const gmx_simd_float_t VD1 = gmx_simd_set1_f(0.43336185284710920150f);
1282 const gmx_simd_float_t VD0 = gmx_simd_set1_f(1.0f);
1284 gmx_simd_float_t z4;
1285 gmx_simd_float_t polyVN0, polyVN1, polyVD0, polyVD1;
1287 z4 = gmx_simd_mul_f(z2, z2);
1289 polyVD1 = gmx_simd_fmadd_f(VD3, z4, VD1);
1290 polyVD0 = gmx_simd_fmadd_f(VD2, z4, VD0);
1291 polyVD0 = gmx_simd_fmadd_f(polyVD1, z2, polyVD0);
1293 polyVD0 = gmx_simd_inv_f(polyVD0);
1295 polyVN0 = gmx_simd_fmadd_f(VN6, z4, VN4);
1296 polyVN1 = gmx_simd_fmadd_f(VN5, z4, VN3);
1297 polyVN0 = gmx_simd_fmadd_f(polyVN0, z4, VN2);
1298 polyVN1 = gmx_simd_fmadd_f(polyVN1, z4, VN1);
1299 polyVN0 = gmx_simd_fmadd_f(polyVN0, z4, VN0);
1300 polyVN0 = gmx_simd_fmadd_f(polyVN1, z2, polyVN0);
1302 return gmx_simd_mul_f(polyVN0, polyVD0);
1308 #ifdef GMX_SIMD_HAVE_DOUBLE
1310 /*! \name Double precision SIMD math functions
1312 * \note In most cases you should use the real-precision functions instead.
1316 /****************************************
1317 * DOUBLE PRECISION SIMD MATH FUNCTIONS *
1318 ****************************************/
1320 /*! \brief SIMD utility function to sum a+b+c+d for SIMD doubles.
1322 * \copydetails gmx_simd_sum4_f
1324 static gmx_inline gmx_simd_double_t gmx_simdcall
1325 gmx_simd_sum4_d(gmx_simd_double_t a, gmx_simd_double_t b,
1326 gmx_simd_double_t c, gmx_simd_double_t d)
1328 return gmx_simd_add_d(gmx_simd_add_d(a, b), gmx_simd_add_d(c, d));
1331 /*! \brief Return -a if b is negative, SIMD double.
1333 * You should normally call the real-precision routine \ref gmx_simd_xor_sign_r.
1335 * \param a Values to set sign for
1336 * \param b Values used to set sign
1337 * \return if b is negative, the sign of a will be changed.
1339 * This is equivalent to doing an xor operation on a with the sign bit of b,
1340 * with the exception that negative zero is not considered to be negative
1341 * on architectures where \ref GMX_SIMD_HAVE_LOGICAL is not set.
1343 static gmx_inline gmx_simd_double_t gmx_simdcall
1344 gmx_simd_xor_sign_d(gmx_simd_double_t a, gmx_simd_double_t b)
1346 #ifdef GMX_SIMD_HAVE_LOGICAL
1347 return gmx_simd_xor_d(a, gmx_simd_and_d(gmx_simd_set1_d(GMX_DOUBLE_NEGZERO), b));
1349 return gmx_simd_blendv_d(a, gmx_simd_fneg_d(a), gmx_simd_cmplt_d(b, gmx_simd_setzero_d()));
1353 /*! \brief Perform one Newton-Raphson iteration to improve 1/sqrt(x) for SIMD double.
1355 * \copydetails gmx_simd_rsqrt_iter_f
1357 static gmx_inline gmx_simd_double_t gmx_simdcall
1358 gmx_simd_rsqrt_iter_d(gmx_simd_double_t lu, gmx_simd_double_t x)
1360 #ifdef GMX_SIMD_HAVE_FMA
1361 return gmx_simd_fmadd_d(gmx_simd_fnmadd_d(x, gmx_simd_mul_d(lu, lu), gmx_simd_set1_d(1.0)), gmx_simd_mul_d(lu, gmx_simd_set1_d(0.5)), lu);
1363 return gmx_simd_mul_d(gmx_simd_set1_d(0.5), gmx_simd_mul_d(gmx_simd_sub_d(gmx_simd_set1_d(3.0), gmx_simd_mul_d(gmx_simd_mul_d(lu, lu), x)), lu));
1368 /*! \brief Calculate 1/sqrt(x) for SIMD double
1370 * \copydetails gmx_simd_invsqrt_f
1372 static gmx_inline gmx_simd_double_t gmx_simdcall
1373 gmx_simd_invsqrt_d(gmx_simd_double_t x)
1375 gmx_simd_double_t lu = gmx_simd_rsqrt_d(x);
1376 #if (GMX_SIMD_RSQRT_BITS < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
1377 lu = gmx_simd_rsqrt_iter_d(lu, x);
1379 #if (GMX_SIMD_RSQRT_BITS*2 < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
1380 lu = gmx_simd_rsqrt_iter_d(lu, x);
1382 #if (GMX_SIMD_RSQRT_BITS*4 < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
1383 lu = gmx_simd_rsqrt_iter_d(lu, x);
1385 #if (GMX_SIMD_RSQRT_BITS*8 < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
1386 lu = gmx_simd_rsqrt_iter_d(lu, x);
1391 /*! \brief Calculate 1/sqrt(x) for two SIMD doubles.
1393 * \copydetails gmx_simd_invsqrt_pair_f
1395 static gmx_inline void gmx_simdcall
1396 gmx_simd_invsqrt_pair_d(gmx_simd_double_t x0, gmx_simd_double_t x1,
1397 gmx_simd_double_t *out0, gmx_simd_double_t *out1)
1399 #if (defined GMX_SIMD_HAVE_FLOAT) && (GMX_SIMD_FLOAT_WIDTH == 2*GMX_SIMD_DOUBLE_WIDTH) && (GMX_SIMD_RSQRT_BITS < 22)
1400 gmx_simd_float_t xf = gmx_simd_cvt_dd2f(x0, x1);
1401 gmx_simd_float_t luf = gmx_simd_rsqrt_f(xf);
1402 gmx_simd_double_t lu0, lu1;
1403 /* Intermediate target is single - mantissa+1 bits */
1404 #if (GMX_SIMD_RSQRT_BITS < GMX_SIMD_MATH_TARGET_SINGLE_BITS)
1405 luf = gmx_simd_rsqrt_iter_f(luf, xf);
1407 #if (GMX_SIMD_RSQRT_BITS*2 < GMX_SIMD_MATH_TARGET_SINGLE_BITS)
1408 luf = gmx_simd_rsqrt_iter_f(luf, xf);
1410 #if (GMX_SIMD_RSQRT_BITS*4 < GMX_SIMD_MATH_TARGET_SINGLE_BITS)
1411 luf = gmx_simd_rsqrt_iter_f(luf, xf);
1413 gmx_simd_cvt_f2dd(luf, &lu0, &lu1);
1414 /* Last iteration(s) performed in double - if we had 22 bits, this gets us to 44 (~1e-15) */
1415 #if (GMX_SIMD_MATH_TARGET_SINGLE_BITS < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
1416 lu0 = gmx_simd_rsqrt_iter_d(lu0, x0);
1417 lu1 = gmx_simd_rsqrt_iter_d(lu1, x1);
1419 #if (GMX_SIMD_MATH_TARGET_SINGLE_BITS*2 < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
1420 lu0 = gmx_simd_rsqrt_iter_d(lu0, x0);
1421 lu1 = gmx_simd_rsqrt_iter_d(lu1, x1);
1426 *out0 = gmx_simd_invsqrt_d(x0);
1427 *out1 = gmx_simd_invsqrt_d(x1);
1431 /*! \brief Perform one Newton-Raphson iteration to improve 1/x for SIMD double.
1433 * \copydetails gmx_simd_rcp_iter_f
1435 static gmx_inline gmx_simd_double_t gmx_simdcall
1436 gmx_simd_rcp_iter_d(gmx_simd_double_t lu, gmx_simd_double_t x)
1438 return gmx_simd_mul_d(lu, gmx_simd_fnmadd_d(lu, x, gmx_simd_set1_d(2.0)));
1441 /*! \brief Calculate 1/x for SIMD double.
1443 * \copydetails gmx_simd_inv_f
1445 static gmx_inline gmx_simd_double_t gmx_simdcall
1446 gmx_simd_inv_d(gmx_simd_double_t x)
1448 gmx_simd_double_t lu = gmx_simd_rcp_d(x);
1449 #if (GMX_SIMD_RCP_BITS < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
1450 lu = gmx_simd_rcp_iter_d(lu, x);
1452 #if (GMX_SIMD_RCP_BITS*2 < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
1453 lu = gmx_simd_rcp_iter_d(lu, x);
1455 #if (GMX_SIMD_RCP_BITS*4 < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
1456 lu = gmx_simd_rcp_iter_d(lu, x);
1458 #if (GMX_SIMD_RCP_BITS*8 < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
1459 lu = gmx_simd_rcp_iter_d(lu, x);
1464 /*! \brief Calculate sqrt(x) correctly for SIMD doubles, including argument 0.0.
1466 * \copydetails gmx_simd_sqrt_f
1468 static gmx_inline gmx_simd_double_t gmx_simdcall
1469 gmx_simd_sqrt_d(gmx_simd_double_t x)
1471 gmx_simd_dbool_t mask;
1472 gmx_simd_double_t res;
1474 mask = gmx_simd_cmpeq_d(x, gmx_simd_setzero_d());
1475 res = gmx_simd_blendnotzero_d(gmx_simd_invsqrt_d(x), mask);
1476 return gmx_simd_mul_d(res, x);
1479 /*! \brief SIMD double log(x). This is the natural logarithm.
1481 * \copydetails gmx_simd_log_f
1483 static gmx_inline gmx_simd_double_t gmx_simdcall
1484 gmx_simd_log_d(gmx_simd_double_t x)
1486 const gmx_simd_double_t half = gmx_simd_set1_d(0.5);
1487 const gmx_simd_double_t one = gmx_simd_set1_d(1.0);
1488 const gmx_simd_double_t sqrt2 = gmx_simd_set1_d(sqrt(2.0));
1489 const gmx_simd_double_t corr = gmx_simd_set1_d(0.693147180559945286226764);
1490 const gmx_simd_double_t CL15 = gmx_simd_set1_d(0.148197055177935105296783);
1491 const gmx_simd_double_t CL13 = gmx_simd_set1_d(0.153108178020442575739679);
1492 const gmx_simd_double_t CL11 = gmx_simd_set1_d(0.181837339521549679055568);
1493 const gmx_simd_double_t CL9 = gmx_simd_set1_d(0.22222194152736701733275);
1494 const gmx_simd_double_t CL7 = gmx_simd_set1_d(0.285714288030134544449368);
1495 const gmx_simd_double_t CL5 = gmx_simd_set1_d(0.399999999989941956712869);
1496 const gmx_simd_double_t CL3 = gmx_simd_set1_d(0.666666666666685503450651);
1497 const gmx_simd_double_t CL1 = gmx_simd_set1_d(2.0);
1498 gmx_simd_double_t fexp, x2, p;
1499 gmx_simd_dbool_t mask;
1501 fexp = gmx_simd_get_exponent_d(x);
1502 x = gmx_simd_get_mantissa_d(x);
1504 mask = gmx_simd_cmplt_d(sqrt2, x);
1505 /* Adjust to non-IEEE format for x>sqrt(2): exponent += 1, mantissa *= 0.5 */
1506 fexp = gmx_simd_add_d(fexp, gmx_simd_blendzero_d(one, mask));
1507 x = gmx_simd_mul_d(x, gmx_simd_blendv_d(one, half, mask));
1509 x = gmx_simd_mul_d( gmx_simd_sub_d(x, one), gmx_simd_inv_d( gmx_simd_add_d(x, one) ) );
1510 x2 = gmx_simd_mul_d(x, x);
1512 p = gmx_simd_fmadd_d(CL15, x2, CL13);
1513 p = gmx_simd_fmadd_d(p, x2, CL11);
1514 p = gmx_simd_fmadd_d(p, x2, CL9);
1515 p = gmx_simd_fmadd_d(p, x2, CL7);
1516 p = gmx_simd_fmadd_d(p, x2, CL5);
1517 p = gmx_simd_fmadd_d(p, x2, CL3);
1518 p = gmx_simd_fmadd_d(p, x2, CL1);
1519 p = gmx_simd_fmadd_d(p, x, gmx_simd_mul_d(corr, fexp));
1524 /*! \brief SIMD double 2^x.
1526 * \copydetails gmx_simd_exp2_f
1528 static gmx_inline gmx_simd_double_t gmx_simdcall
1529 gmx_simd_exp2_d(gmx_simd_double_t x)
1531 const gmx_simd_double_t arglimit = gmx_simd_set1_d(1022.0);
1532 const gmx_simd_double_t CE11 = gmx_simd_set1_d(4.435280790452730022081181e-10);
1533 const gmx_simd_double_t CE10 = gmx_simd_set1_d(7.074105630863314448024247e-09);
1534 const gmx_simd_double_t CE9 = gmx_simd_set1_d(1.017819803432096698472621e-07);
1535 const gmx_simd_double_t CE8 = gmx_simd_set1_d(1.321543308956718799557863e-06);
1536 const gmx_simd_double_t CE7 = gmx_simd_set1_d(0.00001525273348995851746990884);
1537 const gmx_simd_double_t CE6 = gmx_simd_set1_d(0.0001540353046251466849082632);
1538 const gmx_simd_double_t CE5 = gmx_simd_set1_d(0.001333355814678995257307880);
1539 const gmx_simd_double_t CE4 = gmx_simd_set1_d(0.009618129107588335039176502);
1540 const gmx_simd_double_t CE3 = gmx_simd_set1_d(0.05550410866481992147457793);
1541 const gmx_simd_double_t CE2 = gmx_simd_set1_d(0.2402265069591015620470894);
1542 const gmx_simd_double_t CE1 = gmx_simd_set1_d(0.6931471805599453304615075);
1543 const gmx_simd_double_t one = gmx_simd_set1_d(1.0);
1544 gmx_simd_double_t fexppart;
1545 gmx_simd_double_t intpart;
1546 gmx_simd_double_t p;
1547 gmx_simd_dbool_t valuemask;
1549 fexppart = gmx_simd_set_exponent_d(x); /* rounds to nearest int internally */
1550 intpart = gmx_simd_round_d(x); /* use same rounding mode here */
1551 valuemask = gmx_simd_cmple_d(gmx_simd_fabs_d(x), arglimit);
1552 fexppart = gmx_simd_blendzero_d(fexppart, valuemask);
1553 x = gmx_simd_sub_d(x, intpart);
1555 p = gmx_simd_fmadd_d(CE11, x, CE10);
1556 p = gmx_simd_fmadd_d(p, x, CE9);
1557 p = gmx_simd_fmadd_d(p, x, CE8);
1558 p = gmx_simd_fmadd_d(p, x, CE7);
1559 p = gmx_simd_fmadd_d(p, x, CE6);
1560 p = gmx_simd_fmadd_d(p, x, CE5);
1561 p = gmx_simd_fmadd_d(p, x, CE4);
1562 p = gmx_simd_fmadd_d(p, x, CE3);
1563 p = gmx_simd_fmadd_d(p, x, CE2);
1564 p = gmx_simd_fmadd_d(p, x, CE1);
1565 p = gmx_simd_fmadd_d(p, x, one);
1566 x = gmx_simd_mul_d(p, fexppart);
1570 /*! \brief SIMD double exp(x).
1572 * \copydetails gmx_simd_exp_f
1574 static gmx_inline gmx_simd_double_t gmx_simdcall
1575 gmx_simd_exp_d(gmx_simd_double_t x)
1577 const gmx_simd_double_t argscale = gmx_simd_set1_d(1.44269504088896340735992468100);
1578 const gmx_simd_double_t arglimit = gmx_simd_set1_d(1022.0);
1579 const gmx_simd_double_t invargscale0 = gmx_simd_set1_d(0.69314718055966295651160180568695068359375);
1580 const gmx_simd_double_t invargscale1 = gmx_simd_set1_d(2.8235290563031577122588448175013436025525412068e-13);
1581 const gmx_simd_double_t CE12 = gmx_simd_set1_d(2.078375306791423699350304e-09);
1582 const gmx_simd_double_t CE11 = gmx_simd_set1_d(2.518173854179933105218635e-08);
1583 const gmx_simd_double_t CE10 = gmx_simd_set1_d(2.755842049600488770111608e-07);
1584 const gmx_simd_double_t CE9 = gmx_simd_set1_d(2.755691815216689746619849e-06);
1585 const gmx_simd_double_t CE8 = gmx_simd_set1_d(2.480158383706245033920920e-05);
1586 const gmx_simd_double_t CE7 = gmx_simd_set1_d(0.0001984127043518048611841321);
1587 const gmx_simd_double_t CE6 = gmx_simd_set1_d(0.001388888889360258341755930);
1588 const gmx_simd_double_t CE5 = gmx_simd_set1_d(0.008333333332907368102819109);
1589 const gmx_simd_double_t CE4 = gmx_simd_set1_d(0.04166666666663836745814631);
1590 const gmx_simd_double_t CE3 = gmx_simd_set1_d(0.1666666666666796929434570);
1591 const gmx_simd_double_t CE2 = gmx_simd_set1_d(0.5);
1592 const gmx_simd_double_t one = gmx_simd_set1_d(1.0);
1593 gmx_simd_double_t fexppart;
1594 gmx_simd_double_t intpart;
1595 gmx_simd_double_t y, p;
1596 gmx_simd_dbool_t valuemask;
1598 y = gmx_simd_mul_d(x, argscale);
1599 fexppart = gmx_simd_set_exponent_d(y); /* rounds to nearest int internally */
1600 intpart = gmx_simd_round_d(y); /* use same rounding mode here */
1601 valuemask = gmx_simd_cmple_d(gmx_simd_fabs_d(y), arglimit);
1602 fexppart = gmx_simd_blendzero_d(fexppart, valuemask);
1604 /* Extended precision arithmetics */
1605 x = gmx_simd_fnmadd_d(invargscale0, intpart, x);
1606 x = gmx_simd_fnmadd_d(invargscale1, intpart, x);
1608 p = gmx_simd_fmadd_d(CE12, x, CE11);
1609 p = gmx_simd_fmadd_d(p, x, CE10);
1610 p = gmx_simd_fmadd_d(p, x, CE9);
1611 p = gmx_simd_fmadd_d(p, x, CE8);
1612 p = gmx_simd_fmadd_d(p, x, CE7);
1613 p = gmx_simd_fmadd_d(p, x, CE6);
1614 p = gmx_simd_fmadd_d(p, x, CE5);
1615 p = gmx_simd_fmadd_d(p, x, CE4);
1616 p = gmx_simd_fmadd_d(p, x, CE3);
1617 p = gmx_simd_fmadd_d(p, x, CE2);
1618 p = gmx_simd_fmadd_d(p, gmx_simd_mul_d(x, x), gmx_simd_add_d(x, one));
1619 x = gmx_simd_mul_d(p, fexppart);
1623 /*! \brief SIMD double erf(x).
1625 * \copydetails gmx_simd_erf_f
1627 static gmx_inline gmx_simd_double_t gmx_simdcall
1628 gmx_simd_erf_d(gmx_simd_double_t x)
1630 /* Coefficients for minimax approximation of erf(x)=x*(CAoffset + P(x^2)/Q(x^2)) in range [-0.75,0.75] */
1631 const gmx_simd_double_t CAP4 = gmx_simd_set1_d(-0.431780540597889301512e-4);
1632 const gmx_simd_double_t CAP3 = gmx_simd_set1_d(-0.00578562306260059236059);
1633 const gmx_simd_double_t CAP2 = gmx_simd_set1_d(-0.028593586920219752446);
1634 const gmx_simd_double_t CAP1 = gmx_simd_set1_d(-0.315924962948621698209);
1635 const gmx_simd_double_t CAP0 = gmx_simd_set1_d(0.14952975608477029151);
1637 const gmx_simd_double_t CAQ5 = gmx_simd_set1_d(-0.374089300177174709737e-5);
1638 const gmx_simd_double_t CAQ4 = gmx_simd_set1_d(0.00015126584532155383535);
1639 const gmx_simd_double_t CAQ3 = gmx_simd_set1_d(0.00536692680669480725423);
1640 const gmx_simd_double_t CAQ2 = gmx_simd_set1_d(0.0668686825594046122636);
1641 const gmx_simd_double_t CAQ1 = gmx_simd_set1_d(0.402604990869284362773);
1643 const gmx_simd_double_t CAoffset = gmx_simd_set1_d(0.9788494110107421875);
1645 /* Coefficients for minimax approximation of erfc(x)=exp(-x^2)*x*(P(x-1)/Q(x-1)) in range [1.0,4.5] */
1646 const gmx_simd_double_t CBP6 = gmx_simd_set1_d(2.49650423685462752497647637088e-10);
1647 const gmx_simd_double_t CBP5 = gmx_simd_set1_d(0.00119770193298159629350136085658);
1648 const gmx_simd_double_t CBP4 = gmx_simd_set1_d(0.0164944422378370965881008942733);
1649 const gmx_simd_double_t CBP3 = gmx_simd_set1_d(0.0984581468691775932063932439252);
1650 const gmx_simd_double_t CBP2 = gmx_simd_set1_d(0.317364595806937763843589437418);
1651 const gmx_simd_double_t CBP1 = gmx_simd_set1_d(0.554167062641455850932670067075);
1652 const gmx_simd_double_t CBP0 = gmx_simd_set1_d(0.427583576155807163756925301060);
1653 const gmx_simd_double_t CBQ7 = gmx_simd_set1_d(0.00212288829699830145976198384930);
1654 const gmx_simd_double_t CBQ6 = gmx_simd_set1_d(0.0334810979522685300554606393425);
1655 const gmx_simd_double_t CBQ5 = gmx_simd_set1_d(0.2361713785181450957579508850717);
1656 const gmx_simd_double_t CBQ4 = gmx_simd_set1_d(0.955364736493055670530981883072);
1657 const gmx_simd_double_t CBQ3 = gmx_simd_set1_d(2.36815675631420037315349279199);
1658 const gmx_simd_double_t CBQ2 = gmx_simd_set1_d(3.55261649184083035537184223542);
1659 const gmx_simd_double_t CBQ1 = gmx_simd_set1_d(2.93501136050160872574376997993);
1662 /* Coefficients for minimax approximation of erfc(x)=exp(-x^2)/x*(P(1/x)/Q(1/x)) in range [4.5,inf] */
1663 const gmx_simd_double_t CCP6 = gmx_simd_set1_d(-2.8175401114513378771);
1664 const gmx_simd_double_t CCP5 = gmx_simd_set1_d(-3.22729451764143718517);
1665 const gmx_simd_double_t CCP4 = gmx_simd_set1_d(-2.5518551727311523996);
1666 const gmx_simd_double_t CCP3 = gmx_simd_set1_d(-0.687717681153649930619);
1667 const gmx_simd_double_t CCP2 = gmx_simd_set1_d(-0.212652252872804219852);
1668 const gmx_simd_double_t CCP1 = gmx_simd_set1_d(0.0175389834052493308818);
1669 const gmx_simd_double_t CCP0 = gmx_simd_set1_d(0.00628057170626964891937);
1671 const gmx_simd_double_t CCQ6 = gmx_simd_set1_d(5.48409182238641741584);
1672 const gmx_simd_double_t CCQ5 = gmx_simd_set1_d(13.5064170191802889145);
1673 const gmx_simd_double_t CCQ4 = gmx_simd_set1_d(22.9367376522880577224);
1674 const gmx_simd_double_t CCQ3 = gmx_simd_set1_d(15.930646027911794143);
1675 const gmx_simd_double_t CCQ2 = gmx_simd_set1_d(11.0567237927800161565);
1676 const gmx_simd_double_t CCQ1 = gmx_simd_set1_d(2.79257750980575282228);
1678 const gmx_simd_double_t CCoffset = gmx_simd_set1_d(0.5579090118408203125);
1680 const gmx_simd_double_t one = gmx_simd_set1_d(1.0);
1681 const gmx_simd_double_t two = gmx_simd_set1_d(2.0);
1683 gmx_simd_double_t xabs, x2, x4, t, t2, w, w2;
1684 gmx_simd_double_t PolyAP0, PolyAP1, PolyAQ0, PolyAQ1;
1685 gmx_simd_double_t PolyBP0, PolyBP1, PolyBQ0, PolyBQ1;
1686 gmx_simd_double_t PolyCP0, PolyCP1, PolyCQ0, PolyCQ1;
1687 gmx_simd_double_t res_erf, res_erfcB, res_erfcC, res_erfc, res;
1688 gmx_simd_double_t expmx2;
1689 gmx_simd_dbool_t mask;
1691 /* Calculate erf() */
1692 xabs = gmx_simd_fabs_d(x);
1693 x2 = gmx_simd_mul_d(x, x);
1694 x4 = gmx_simd_mul_d(x2, x2);
1696 PolyAP0 = gmx_simd_mul_d(CAP4, x4);
1697 PolyAP1 = gmx_simd_mul_d(CAP3, x4);
1698 PolyAP0 = gmx_simd_add_d(PolyAP0, CAP2);
1699 PolyAP1 = gmx_simd_add_d(PolyAP1, CAP1);
1700 PolyAP0 = gmx_simd_mul_d(PolyAP0, x4);
1701 PolyAP1 = gmx_simd_mul_d(PolyAP1, x2);
1702 PolyAP0 = gmx_simd_add_d(PolyAP0, CAP0);
1703 PolyAP0 = gmx_simd_add_d(PolyAP0, PolyAP1);
1705 PolyAQ1 = gmx_simd_mul_d(CAQ5, x4);
1706 PolyAQ0 = gmx_simd_mul_d(CAQ4, x4);
1707 PolyAQ1 = gmx_simd_add_d(PolyAQ1, CAQ3);
1708 PolyAQ0 = gmx_simd_add_d(PolyAQ0, CAQ2);
1709 PolyAQ1 = gmx_simd_mul_d(PolyAQ1, x4);
1710 PolyAQ0 = gmx_simd_mul_d(PolyAQ0, x4);
1711 PolyAQ1 = gmx_simd_add_d(PolyAQ1, CAQ1);
1712 PolyAQ0 = gmx_simd_add_d(PolyAQ0, one);
1713 PolyAQ1 = gmx_simd_mul_d(PolyAQ1, x2);
1714 PolyAQ0 = gmx_simd_add_d(PolyAQ0, PolyAQ1);
1716 res_erf = gmx_simd_mul_d(PolyAP0, gmx_simd_inv_d(PolyAQ0));
1717 res_erf = gmx_simd_add_d(CAoffset, res_erf);
1718 res_erf = gmx_simd_mul_d(x, res_erf);
1720 /* Calculate erfc() in range [1,4.5] */
1721 t = gmx_simd_sub_d(xabs, one);
1722 t2 = gmx_simd_mul_d(t, t);
1724 PolyBP0 = gmx_simd_mul_d(CBP6, t2);
1725 PolyBP1 = gmx_simd_mul_d(CBP5, t2);
1726 PolyBP0 = gmx_simd_add_d(PolyBP0, CBP4);
1727 PolyBP1 = gmx_simd_add_d(PolyBP1, CBP3);
1728 PolyBP0 = gmx_simd_mul_d(PolyBP0, t2);
1729 PolyBP1 = gmx_simd_mul_d(PolyBP1, t2);
1730 PolyBP0 = gmx_simd_add_d(PolyBP0, CBP2);
1731 PolyBP1 = gmx_simd_add_d(PolyBP1, CBP1);
1732 PolyBP0 = gmx_simd_mul_d(PolyBP0, t2);
1733 PolyBP1 = gmx_simd_mul_d(PolyBP1, t);
1734 PolyBP0 = gmx_simd_add_d(PolyBP0, CBP0);
1735 PolyBP0 = gmx_simd_add_d(PolyBP0, PolyBP1);
1737 PolyBQ1 = gmx_simd_mul_d(CBQ7, t2);
1738 PolyBQ0 = gmx_simd_mul_d(CBQ6, t2);
1739 PolyBQ1 = gmx_simd_add_d(PolyBQ1, CBQ5);
1740 PolyBQ0 = gmx_simd_add_d(PolyBQ0, CBQ4);
1741 PolyBQ1 = gmx_simd_mul_d(PolyBQ1, t2);
1742 PolyBQ0 = gmx_simd_mul_d(PolyBQ0, t2);
1743 PolyBQ1 = gmx_simd_add_d(PolyBQ1, CBQ3);
1744 PolyBQ0 = gmx_simd_add_d(PolyBQ0, CBQ2);
1745 PolyBQ1 = gmx_simd_mul_d(PolyBQ1, t2);
1746 PolyBQ0 = gmx_simd_mul_d(PolyBQ0, t2);
1747 PolyBQ1 = gmx_simd_add_d(PolyBQ1, CBQ1);
1748 PolyBQ0 = gmx_simd_add_d(PolyBQ0, one);
1749 PolyBQ1 = gmx_simd_mul_d(PolyBQ1, t);
1750 PolyBQ0 = gmx_simd_add_d(PolyBQ0, PolyBQ1);
1752 res_erfcB = gmx_simd_mul_d(PolyBP0, gmx_simd_inv_d(PolyBQ0));
1754 res_erfcB = gmx_simd_mul_d(res_erfcB, xabs);
1756 /* Calculate erfc() in range [4.5,inf] */
1757 w = gmx_simd_inv_d(xabs);
1758 w2 = gmx_simd_mul_d(w, w);
1760 PolyCP0 = gmx_simd_mul_d(CCP6, w2);
1761 PolyCP1 = gmx_simd_mul_d(CCP5, w2);
1762 PolyCP0 = gmx_simd_add_d(PolyCP0, CCP4);
1763 PolyCP1 = gmx_simd_add_d(PolyCP1, CCP3);
1764 PolyCP0 = gmx_simd_mul_d(PolyCP0, w2);
1765 PolyCP1 = gmx_simd_mul_d(PolyCP1, w2);
1766 PolyCP0 = gmx_simd_add_d(PolyCP0, CCP2);
1767 PolyCP1 = gmx_simd_add_d(PolyCP1, CCP1);
1768 PolyCP0 = gmx_simd_mul_d(PolyCP0, w2);
1769 PolyCP1 = gmx_simd_mul_d(PolyCP1, w);
1770 PolyCP0 = gmx_simd_add_d(PolyCP0, CCP0);
1771 PolyCP0 = gmx_simd_add_d(PolyCP0, PolyCP1);
1773 PolyCQ0 = gmx_simd_mul_d(CCQ6, w2);
1774 PolyCQ1 = gmx_simd_mul_d(CCQ5, w2);
1775 PolyCQ0 = gmx_simd_add_d(PolyCQ0, CCQ4);
1776 PolyCQ1 = gmx_simd_add_d(PolyCQ1, CCQ3);
1777 PolyCQ0 = gmx_simd_mul_d(PolyCQ0, w2);
1778 PolyCQ1 = gmx_simd_mul_d(PolyCQ1, w2);
1779 PolyCQ0 = gmx_simd_add_d(PolyCQ0, CCQ2);
1780 PolyCQ1 = gmx_simd_add_d(PolyCQ1, CCQ1);
1781 PolyCQ0 = gmx_simd_mul_d(PolyCQ0, w2);
1782 PolyCQ1 = gmx_simd_mul_d(PolyCQ1, w);
1783 PolyCQ0 = gmx_simd_add_d(PolyCQ0, one);
1784 PolyCQ0 = gmx_simd_add_d(PolyCQ0, PolyCQ1);
1786 expmx2 = gmx_simd_exp_d( gmx_simd_fneg_d(x2) );
1788 res_erfcC = gmx_simd_mul_d(PolyCP0, gmx_simd_inv_d(PolyCQ0));
1789 res_erfcC = gmx_simd_add_d(res_erfcC, CCoffset);
1790 res_erfcC = gmx_simd_mul_d(res_erfcC, w);
1792 mask = gmx_simd_cmplt_d(gmx_simd_set1_d(4.5), xabs);
1793 res_erfc = gmx_simd_blendv_d(res_erfcB, res_erfcC, mask);
1795 res_erfc = gmx_simd_mul_d(res_erfc, expmx2);
1797 /* erfc(x<0) = 2-erfc(|x|) */
1798 mask = gmx_simd_cmplt_d(x, gmx_simd_setzero_d());
1799 res_erfc = gmx_simd_blendv_d(res_erfc, gmx_simd_sub_d(two, res_erfc), mask);
1801 /* Select erf() or erfc() */
1802 mask = gmx_simd_cmplt_d(xabs, one);
1803 res = gmx_simd_blendv_d(gmx_simd_sub_d(one, res_erfc), res_erf, mask);
1808 /*! \brief SIMD double erfc(x).
1810 * \copydetails gmx_simd_erfc_f
1812 static gmx_inline gmx_simd_double_t gmx_simdcall
1813 gmx_simd_erfc_d(gmx_simd_double_t x)
1815 /* Coefficients for minimax approximation of erf(x)=x*(CAoffset + P(x^2)/Q(x^2)) in range [-0.75,0.75] */
1816 const gmx_simd_double_t CAP4 = gmx_simd_set1_d(-0.431780540597889301512e-4);
1817 const gmx_simd_double_t CAP3 = gmx_simd_set1_d(-0.00578562306260059236059);
1818 const gmx_simd_double_t CAP2 = gmx_simd_set1_d(-0.028593586920219752446);
1819 const gmx_simd_double_t CAP1 = gmx_simd_set1_d(-0.315924962948621698209);
1820 const gmx_simd_double_t CAP0 = gmx_simd_set1_d(0.14952975608477029151);
1822 const gmx_simd_double_t CAQ5 = gmx_simd_set1_d(-0.374089300177174709737e-5);
1823 const gmx_simd_double_t CAQ4 = gmx_simd_set1_d(0.00015126584532155383535);
1824 const gmx_simd_double_t CAQ3 = gmx_simd_set1_d(0.00536692680669480725423);
1825 const gmx_simd_double_t CAQ2 = gmx_simd_set1_d(0.0668686825594046122636);
1826 const gmx_simd_double_t CAQ1 = gmx_simd_set1_d(0.402604990869284362773);
1828 const gmx_simd_double_t CAoffset = gmx_simd_set1_d(0.9788494110107421875);
1830 /* Coefficients for minimax approximation of erfc(x)=exp(-x^2)*x*(P(x-1)/Q(x-1)) in range [1.0,4.5] */
1831 const gmx_simd_double_t CBP6 = gmx_simd_set1_d(2.49650423685462752497647637088e-10);
1832 const gmx_simd_double_t CBP5 = gmx_simd_set1_d(0.00119770193298159629350136085658);
1833 const gmx_simd_double_t CBP4 = gmx_simd_set1_d(0.0164944422378370965881008942733);
1834 const gmx_simd_double_t CBP3 = gmx_simd_set1_d(0.0984581468691775932063932439252);
1835 const gmx_simd_double_t CBP2 = gmx_simd_set1_d(0.317364595806937763843589437418);
1836 const gmx_simd_double_t CBP1 = gmx_simd_set1_d(0.554167062641455850932670067075);
1837 const gmx_simd_double_t CBP0 = gmx_simd_set1_d(0.427583576155807163756925301060);
1838 const gmx_simd_double_t CBQ7 = gmx_simd_set1_d(0.00212288829699830145976198384930);
1839 const gmx_simd_double_t CBQ6 = gmx_simd_set1_d(0.0334810979522685300554606393425);
1840 const gmx_simd_double_t CBQ5 = gmx_simd_set1_d(0.2361713785181450957579508850717);
1841 const gmx_simd_double_t CBQ4 = gmx_simd_set1_d(0.955364736493055670530981883072);
1842 const gmx_simd_double_t CBQ3 = gmx_simd_set1_d(2.36815675631420037315349279199);
1843 const gmx_simd_double_t CBQ2 = gmx_simd_set1_d(3.55261649184083035537184223542);
1844 const gmx_simd_double_t CBQ1 = gmx_simd_set1_d(2.93501136050160872574376997993);
1847 /* Coefficients for minimax approximation of erfc(x)=exp(-x^2)/x*(P(1/x)/Q(1/x)) in range [4.5,inf] */
1848 const gmx_simd_double_t CCP6 = gmx_simd_set1_d(-2.8175401114513378771);
1849 const gmx_simd_double_t CCP5 = gmx_simd_set1_d(-3.22729451764143718517);
1850 const gmx_simd_double_t CCP4 = gmx_simd_set1_d(-2.5518551727311523996);
1851 const gmx_simd_double_t CCP3 = gmx_simd_set1_d(-0.687717681153649930619);
1852 const gmx_simd_double_t CCP2 = gmx_simd_set1_d(-0.212652252872804219852);
1853 const gmx_simd_double_t CCP1 = gmx_simd_set1_d(0.0175389834052493308818);
1854 const gmx_simd_double_t CCP0 = gmx_simd_set1_d(0.00628057170626964891937);
1856 const gmx_simd_double_t CCQ6 = gmx_simd_set1_d(5.48409182238641741584);
1857 const gmx_simd_double_t CCQ5 = gmx_simd_set1_d(13.5064170191802889145);
1858 const gmx_simd_double_t CCQ4 = gmx_simd_set1_d(22.9367376522880577224);
1859 const gmx_simd_double_t CCQ3 = gmx_simd_set1_d(15.930646027911794143);
1860 const gmx_simd_double_t CCQ2 = gmx_simd_set1_d(11.0567237927800161565);
1861 const gmx_simd_double_t CCQ1 = gmx_simd_set1_d(2.79257750980575282228);
1863 const gmx_simd_double_t CCoffset = gmx_simd_set1_d(0.5579090118408203125);
1865 const gmx_simd_double_t one = gmx_simd_set1_d(1.0);
1866 const gmx_simd_double_t two = gmx_simd_set1_d(2.0);
1868 gmx_simd_double_t xabs, x2, x4, t, t2, w, w2;
1869 gmx_simd_double_t PolyAP0, PolyAP1, PolyAQ0, PolyAQ1;
1870 gmx_simd_double_t PolyBP0, PolyBP1, PolyBQ0, PolyBQ1;
1871 gmx_simd_double_t PolyCP0, PolyCP1, PolyCQ0, PolyCQ1;
1872 gmx_simd_double_t res_erf, res_erfcB, res_erfcC, res_erfc, res;
1873 gmx_simd_double_t expmx2;
1874 gmx_simd_dbool_t mask;
1876 /* Calculate erf() */
1877 xabs = gmx_simd_fabs_d(x);
1878 x2 = gmx_simd_mul_d(x, x);
1879 x4 = gmx_simd_mul_d(x2, x2);
1881 PolyAP0 = gmx_simd_mul_d(CAP4, x4);
1882 PolyAP1 = gmx_simd_mul_d(CAP3, x4);
1883 PolyAP0 = gmx_simd_add_d(PolyAP0, CAP2);
1884 PolyAP1 = gmx_simd_add_d(PolyAP1, CAP1);
1885 PolyAP0 = gmx_simd_mul_d(PolyAP0, x4);
1886 PolyAP1 = gmx_simd_mul_d(PolyAP1, x2);
1887 PolyAP0 = gmx_simd_add_d(PolyAP0, CAP0);
1888 PolyAP0 = gmx_simd_add_d(PolyAP0, PolyAP1);
1890 PolyAQ1 = gmx_simd_mul_d(CAQ5, x4);
1891 PolyAQ0 = gmx_simd_mul_d(CAQ4, x4);
1892 PolyAQ1 = gmx_simd_add_d(PolyAQ1, CAQ3);
1893 PolyAQ0 = gmx_simd_add_d(PolyAQ0, CAQ2);
1894 PolyAQ1 = gmx_simd_mul_d(PolyAQ1, x4);
1895 PolyAQ0 = gmx_simd_mul_d(PolyAQ0, x4);
1896 PolyAQ1 = gmx_simd_add_d(PolyAQ1, CAQ1);
1897 PolyAQ0 = gmx_simd_add_d(PolyAQ0, one);
1898 PolyAQ1 = gmx_simd_mul_d(PolyAQ1, x2);
1899 PolyAQ0 = gmx_simd_add_d(PolyAQ0, PolyAQ1);
1901 res_erf = gmx_simd_mul_d(PolyAP0, gmx_simd_inv_d(PolyAQ0));
1902 res_erf = gmx_simd_add_d(CAoffset, res_erf);
1903 res_erf = gmx_simd_mul_d(x, res_erf);
1905 /* Calculate erfc() in range [1,4.5] */
1906 t = gmx_simd_sub_d(xabs, one);
1907 t2 = gmx_simd_mul_d(t, t);
1909 PolyBP0 = gmx_simd_mul_d(CBP6, t2);
1910 PolyBP1 = gmx_simd_mul_d(CBP5, t2);
1911 PolyBP0 = gmx_simd_add_d(PolyBP0, CBP4);
1912 PolyBP1 = gmx_simd_add_d(PolyBP1, CBP3);
1913 PolyBP0 = gmx_simd_mul_d(PolyBP0, t2);
1914 PolyBP1 = gmx_simd_mul_d(PolyBP1, t2);
1915 PolyBP0 = gmx_simd_add_d(PolyBP0, CBP2);
1916 PolyBP1 = gmx_simd_add_d(PolyBP1, CBP1);
1917 PolyBP0 = gmx_simd_mul_d(PolyBP0, t2);
1918 PolyBP1 = gmx_simd_mul_d(PolyBP1, t);
1919 PolyBP0 = gmx_simd_add_d(PolyBP0, CBP0);
1920 PolyBP0 = gmx_simd_add_d(PolyBP0, PolyBP1);
1922 PolyBQ1 = gmx_simd_mul_d(CBQ7, t2);
1923 PolyBQ0 = gmx_simd_mul_d(CBQ6, t2);
1924 PolyBQ1 = gmx_simd_add_d(PolyBQ1, CBQ5);
1925 PolyBQ0 = gmx_simd_add_d(PolyBQ0, CBQ4);
1926 PolyBQ1 = gmx_simd_mul_d(PolyBQ1, t2);
1927 PolyBQ0 = gmx_simd_mul_d(PolyBQ0, t2);
1928 PolyBQ1 = gmx_simd_add_d(PolyBQ1, CBQ3);
1929 PolyBQ0 = gmx_simd_add_d(PolyBQ0, CBQ2);
1930 PolyBQ1 = gmx_simd_mul_d(PolyBQ1, t2);
1931 PolyBQ0 = gmx_simd_mul_d(PolyBQ0, t2);
1932 PolyBQ1 = gmx_simd_add_d(PolyBQ1, CBQ1);
1933 PolyBQ0 = gmx_simd_add_d(PolyBQ0, one);
1934 PolyBQ1 = gmx_simd_mul_d(PolyBQ1, t);
1935 PolyBQ0 = gmx_simd_add_d(PolyBQ0, PolyBQ1);
1937 res_erfcB = gmx_simd_mul_d(PolyBP0, gmx_simd_inv_d(PolyBQ0));
1939 res_erfcB = gmx_simd_mul_d(res_erfcB, xabs);
1941 /* Calculate erfc() in range [4.5,inf] */
1942 w = gmx_simd_inv_d(xabs);
1943 w2 = gmx_simd_mul_d(w, w);
1945 PolyCP0 = gmx_simd_mul_d(CCP6, w2);
1946 PolyCP1 = gmx_simd_mul_d(CCP5, w2);
1947 PolyCP0 = gmx_simd_add_d(PolyCP0, CCP4);
1948 PolyCP1 = gmx_simd_add_d(PolyCP1, CCP3);
1949 PolyCP0 = gmx_simd_mul_d(PolyCP0, w2);
1950 PolyCP1 = gmx_simd_mul_d(PolyCP1, w2);
1951 PolyCP0 = gmx_simd_add_d(PolyCP0, CCP2);
1952 PolyCP1 = gmx_simd_add_d(PolyCP1, CCP1);
1953 PolyCP0 = gmx_simd_mul_d(PolyCP0, w2);
1954 PolyCP1 = gmx_simd_mul_d(PolyCP1, w);
1955 PolyCP0 = gmx_simd_add_d(PolyCP0, CCP0);
1956 PolyCP0 = gmx_simd_add_d(PolyCP0, PolyCP1);
1958 PolyCQ0 = gmx_simd_mul_d(CCQ6, w2);
1959 PolyCQ1 = gmx_simd_mul_d(CCQ5, w2);
1960 PolyCQ0 = gmx_simd_add_d(PolyCQ0, CCQ4);
1961 PolyCQ1 = gmx_simd_add_d(PolyCQ1, CCQ3);
1962 PolyCQ0 = gmx_simd_mul_d(PolyCQ0, w2);
1963 PolyCQ1 = gmx_simd_mul_d(PolyCQ1, w2);
1964 PolyCQ0 = gmx_simd_add_d(PolyCQ0, CCQ2);
1965 PolyCQ1 = gmx_simd_add_d(PolyCQ1, CCQ1);
1966 PolyCQ0 = gmx_simd_mul_d(PolyCQ0, w2);
1967 PolyCQ1 = gmx_simd_mul_d(PolyCQ1, w);
1968 PolyCQ0 = gmx_simd_add_d(PolyCQ0, one);
1969 PolyCQ0 = gmx_simd_add_d(PolyCQ0, PolyCQ1);
1971 expmx2 = gmx_simd_exp_d( gmx_simd_fneg_d(x2) );
1973 res_erfcC = gmx_simd_mul_d(PolyCP0, gmx_simd_inv_d(PolyCQ0));
1974 res_erfcC = gmx_simd_add_d(res_erfcC, CCoffset);
1975 res_erfcC = gmx_simd_mul_d(res_erfcC, w);
1977 mask = gmx_simd_cmplt_d(gmx_simd_set1_d(4.5), xabs);
1978 res_erfc = gmx_simd_blendv_d(res_erfcB, res_erfcC, mask);
1980 res_erfc = gmx_simd_mul_d(res_erfc, expmx2);
1982 /* erfc(x<0) = 2-erfc(|x|) */
1983 mask = gmx_simd_cmplt_d(x, gmx_simd_setzero_d());
1984 res_erfc = gmx_simd_blendv_d(res_erfc, gmx_simd_sub_d(two, res_erfc), mask);
1986 /* Select erf() or erfc() */
1987 mask = gmx_simd_cmplt_d(xabs, one);
1988 res = gmx_simd_blendv_d(res_erfc, gmx_simd_sub_d(one, res_erf), mask);
1993 /*! \brief SIMD double sin \& cos.
1995 * \copydetails gmx_simd_sincos_f
1997 static gmx_inline void gmx_simdcall
1998 gmx_simd_sincos_d(gmx_simd_double_t x, gmx_simd_double_t *sinval, gmx_simd_double_t *cosval)
2000 /* Constants to subtract Pi/4*x from y while minimizing precision loss */
2001 const gmx_simd_double_t argred0 = gmx_simd_set1_d(2*0.78539816290140151978);
2002 const gmx_simd_double_t argred1 = gmx_simd_set1_d(2*4.9604678871439933374e-10);
2003 const gmx_simd_double_t argred2 = gmx_simd_set1_d(2*1.1258708853173288931e-18);
2004 const gmx_simd_double_t argred3 = gmx_simd_set1_d(2*1.7607799325916000908e-27);
2005 const gmx_simd_double_t two_over_pi = gmx_simd_set1_d(2.0/M_PI);
2006 const gmx_simd_double_t const_sin5 = gmx_simd_set1_d( 1.58938307283228937328511e-10);
2007 const gmx_simd_double_t const_sin4 = gmx_simd_set1_d(-2.50506943502539773349318e-08);
2008 const gmx_simd_double_t const_sin3 = gmx_simd_set1_d( 2.75573131776846360512547e-06);
2009 const gmx_simd_double_t const_sin2 = gmx_simd_set1_d(-0.000198412698278911770864914);
2010 const gmx_simd_double_t const_sin1 = gmx_simd_set1_d( 0.0083333333333191845961746);
2011 const gmx_simd_double_t const_sin0 = gmx_simd_set1_d(-0.166666666666666130709393);
2013 const gmx_simd_double_t const_cos7 = gmx_simd_set1_d(-1.13615350239097429531523e-11);
2014 const gmx_simd_double_t const_cos6 = gmx_simd_set1_d( 2.08757471207040055479366e-09);
2015 const gmx_simd_double_t const_cos5 = gmx_simd_set1_d(-2.75573144028847567498567e-07);
2016 const gmx_simd_double_t const_cos4 = gmx_simd_set1_d( 2.48015872890001867311915e-05);
2017 const gmx_simd_double_t const_cos3 = gmx_simd_set1_d(-0.00138888888888714019282329);
2018 const gmx_simd_double_t const_cos2 = gmx_simd_set1_d( 0.0416666666666665519592062);
2019 const gmx_simd_double_t half = gmx_simd_set1_d(0.5);
2020 const gmx_simd_double_t one = gmx_simd_set1_d(1.0);
2021 gmx_simd_double_t ssign, csign;
2022 gmx_simd_double_t x2, y, z, psin, pcos, sss, ccc;
2023 gmx_simd_dbool_t mask;
2024 #if (defined GMX_SIMD_HAVE_DINT32) && (defined GMX_SIMD_HAVE_DINT32_ARITHMETICS) && (defined GMX_SIMD_HAVE_LOGICAL)
2025 const gmx_simd_dint32_t ione = gmx_simd_set1_di(1);
2026 const gmx_simd_dint32_t itwo = gmx_simd_set1_di(2);
2027 gmx_simd_dint32_t iy;
2029 z = gmx_simd_mul_d(x, two_over_pi);
2030 iy = gmx_simd_cvt_d2i(z);
2031 y = gmx_simd_round_d(z);
2033 mask = gmx_simd_cvt_dib2db(gmx_simd_cmpeq_di(gmx_simd_and_di(iy, ione), gmx_simd_setzero_di()));
2034 ssign = gmx_simd_blendzero_d(gmx_simd_set1_d(GMX_DOUBLE_NEGZERO), gmx_simd_cvt_dib2db(gmx_simd_cmpeq_di(gmx_simd_and_di(iy, itwo), itwo)));
2035 csign = gmx_simd_blendzero_d(gmx_simd_set1_d(GMX_DOUBLE_NEGZERO), gmx_simd_cvt_dib2db(gmx_simd_cmpeq_di(gmx_simd_and_di(gmx_simd_add_di(iy, ione), itwo), itwo)));
2037 const gmx_simd_double_t quarter = gmx_simd_set1_d(0.25);
2038 const gmx_simd_double_t minusquarter = gmx_simd_set1_d(-0.25);
2039 gmx_simd_double_t q;
2040 gmx_simd_dbool_t m1, m2, m3;
2042 /* The most obvious way to find the arguments quadrant in the unit circle
2043 * to calculate the sign is to use integer arithmetic, but that is not
2044 * present in all SIMD implementations. As an alternative, we have devised a
2045 * pure floating-point algorithm that uses truncation for argument reduction
2046 * so that we get a new value 0<=q<1 over the unit circle, and then
2047 * do floating-point comparisons with fractions. This is likely to be
2048 * slightly slower (~10%) due to the longer latencies of floating-point, so
2049 * we only use it when integer SIMD arithmetic is not present.
2052 x = gmx_simd_fabs_d(x);
2053 /* It is critical that half-way cases are rounded down */
2054 z = gmx_simd_fmadd_d(x, two_over_pi, half);
2055 y = gmx_simd_trunc_d(z);
2056 q = gmx_simd_mul_d(z, quarter);
2057 q = gmx_simd_sub_d(q, gmx_simd_trunc_d(q));
2058 /* z now starts at 0.0 for x=-pi/4 (although neg. values cannot occur), and
2059 * then increased by 1.0 as x increases by 2*Pi, when it resets to 0.0.
2060 * This removes the 2*Pi periodicity without using any integer arithmetic.
2061 * First check if y had the value 2 or 3, set csign if true.
2063 q = gmx_simd_sub_d(q, half);
2064 /* If we have logical operations we can work directly on the signbit, which
2065 * saves instructions. Otherwise we need to represent signs as +1.0/-1.0.
2066 * Thus, if you are altering defines to debug alternative code paths, the
2067 * two GMX_SIMD_HAVE_LOGICAL sections in this routine must either both be
2068 * active or inactive - you will get errors if only one is used.
2070 # ifdef GMX_SIMD_HAVE_LOGICAL
2071 ssign = gmx_simd_and_d(ssign, gmx_simd_set1_d(GMX_DOUBLE_NEGZERO));
2072 csign = gmx_simd_andnot_d(q, gmx_simd_set1_d(GMX_DOUBLE_NEGZERO));
2073 ssign = gmx_simd_xor_d(ssign, csign);
2075 csign = gmx_simd_xor_sign_d(gmx_simd_set1_d(-1.0), q);
2076 ssign = gmx_simd_xor_sign_d(ssign, csign); /* swap ssign if csign was set. */
2078 /* Check if y had value 1 or 3 (remember we subtracted 0.5 from q) */
2079 m1 = gmx_simd_cmplt_d(q, minusquarter);
2080 m2 = gmx_simd_cmple_d(gmx_simd_setzero_d(), q);
2081 m3 = gmx_simd_cmplt_d(q, quarter);
2082 m2 = gmx_simd_and_db(m2, m3);
2083 mask = gmx_simd_or_db(m1, m2);
2084 /* where mask is FALSE, set sign. */
2085 csign = gmx_simd_xor_sign_d(csign, gmx_simd_blendv_d(gmx_simd_set1_d(-1.0), one, mask));
2087 x = gmx_simd_fnmadd_d(y, argred0, x);
2088 x = gmx_simd_fnmadd_d(y, argred1, x);
2089 x = gmx_simd_fnmadd_d(y, argred2, x);
2090 x = gmx_simd_fnmadd_d(y, argred3, x);
2091 x2 = gmx_simd_mul_d(x, x);
2093 psin = gmx_simd_fmadd_d(const_sin5, x2, const_sin4);
2094 psin = gmx_simd_fmadd_d(psin, x2, const_sin3);
2095 psin = gmx_simd_fmadd_d(psin, x2, const_sin2);
2096 psin = gmx_simd_fmadd_d(psin, x2, const_sin1);
2097 psin = gmx_simd_fmadd_d(psin, x2, const_sin0);
2098 psin = gmx_simd_fmadd_d(psin, gmx_simd_mul_d(x2, x), x);
2100 pcos = gmx_simd_fmadd_d(const_cos7, x2, const_cos6);
2101 pcos = gmx_simd_fmadd_d(pcos, x2, const_cos5);
2102 pcos = gmx_simd_fmadd_d(pcos, x2, const_cos4);
2103 pcos = gmx_simd_fmadd_d(pcos, x2, const_cos3);
2104 pcos = gmx_simd_fmadd_d(pcos, x2, const_cos2);
2105 pcos = gmx_simd_fmsub_d(pcos, x2, half);
2106 pcos = gmx_simd_fmadd_d(pcos, x2, one);
2108 sss = gmx_simd_blendv_d(pcos, psin, mask);
2109 ccc = gmx_simd_blendv_d(psin, pcos, mask);
2110 /* See comment for GMX_SIMD_HAVE_LOGICAL section above. */
2111 #ifdef GMX_SIMD_HAVE_LOGICAL
2112 *sinval = gmx_simd_xor_d(sss, ssign);
2113 *cosval = gmx_simd_xor_d(ccc, csign);
2115 *sinval = gmx_simd_xor_sign_d(sss, ssign);
2116 *cosval = gmx_simd_xor_sign_d(ccc, csign);
2120 /*! \brief SIMD double sin(x).
2122 * \copydetails gmx_simd_sin_f
2124 static gmx_inline gmx_simd_double_t gmx_simdcall
2125 gmx_simd_sin_d(gmx_simd_double_t x)
2127 gmx_simd_double_t s, c;
2128 gmx_simd_sincos_d(x, &s, &c);
2132 /*! \brief SIMD double cos(x).
2134 * \copydetails gmx_simd_cos_f
2136 static gmx_inline gmx_simd_double_t gmx_simdcall
2137 gmx_simd_cos_d(gmx_simd_double_t x)
2139 gmx_simd_double_t s, c;
2140 gmx_simd_sincos_d(x, &s, &c);
2144 /*! \brief SIMD double tan(x).
2146 * \copydetails gmx_simd_tan_f
2148 static gmx_inline gmx_simd_double_t gmx_simdcall
2149 gmx_simd_tan_d(gmx_simd_double_t x)
2151 const gmx_simd_double_t argred0 = gmx_simd_set1_d(2*0.78539816290140151978);
2152 const gmx_simd_double_t argred1 = gmx_simd_set1_d(2*4.9604678871439933374e-10);
2153 const gmx_simd_double_t argred2 = gmx_simd_set1_d(2*1.1258708853173288931e-18);
2154 const gmx_simd_double_t argred3 = gmx_simd_set1_d(2*1.7607799325916000908e-27);
2155 const gmx_simd_double_t two_over_pi = gmx_simd_set1_d(2.0/M_PI);
2156 const gmx_simd_double_t CT15 = gmx_simd_set1_d(1.01419718511083373224408e-05);
2157 const gmx_simd_double_t CT14 = gmx_simd_set1_d(-2.59519791585924697698614e-05);
2158 const gmx_simd_double_t CT13 = gmx_simd_set1_d(5.23388081915899855325186e-05);
2159 const gmx_simd_double_t CT12 = gmx_simd_set1_d(-3.05033014433946488225616e-05);
2160 const gmx_simd_double_t CT11 = gmx_simd_set1_d(7.14707504084242744267497e-05);
2161 const gmx_simd_double_t CT10 = gmx_simd_set1_d(8.09674518280159187045078e-05);
2162 const gmx_simd_double_t CT9 = gmx_simd_set1_d(0.000244884931879331847054404);
2163 const gmx_simd_double_t CT8 = gmx_simd_set1_d(0.000588505168743587154904506);
2164 const gmx_simd_double_t CT7 = gmx_simd_set1_d(0.00145612788922812427978848);
2165 const gmx_simd_double_t CT6 = gmx_simd_set1_d(0.00359208743836906619142924);
2166 const gmx_simd_double_t CT5 = gmx_simd_set1_d(0.00886323944362401618113356);
2167 const gmx_simd_double_t CT4 = gmx_simd_set1_d(0.0218694882853846389592078);
2168 const gmx_simd_double_t CT3 = gmx_simd_set1_d(0.0539682539781298417636002);
2169 const gmx_simd_double_t CT2 = gmx_simd_set1_d(0.133333333333125941821962);
2170 const gmx_simd_double_t CT1 = gmx_simd_set1_d(0.333333333333334980164153);
2172 gmx_simd_double_t x2, p, y, z;
2173 gmx_simd_dbool_t mask;
2175 #if (defined GMX_SIMD_HAVE_DINT32) && (defined GMX_SIMD_HAVE_DINT32_ARITHMETICS) && (defined GMX_SIMD_HAVE_LOGICAL)
2176 gmx_simd_dint32_t iy;
2177 gmx_simd_dint32_t ione = gmx_simd_set1_di(1);
2179 z = gmx_simd_mul_d(x, two_over_pi);
2180 iy = gmx_simd_cvt_d2i(z);
2181 y = gmx_simd_round_d(z);
2182 mask = gmx_simd_cvt_dib2db(gmx_simd_cmpeq_di(gmx_simd_and_di(iy, ione), ione));
2184 x = gmx_simd_fnmadd_d(y, argred0, x);
2185 x = gmx_simd_fnmadd_d(y, argred1, x);
2186 x = gmx_simd_fnmadd_d(y, argred2, x);
2187 x = gmx_simd_fnmadd_d(y, argred3, x);
2188 x = gmx_simd_xor_d(gmx_simd_blendzero_d(gmx_simd_set1_d(GMX_DOUBLE_NEGZERO), mask), x);
2190 const gmx_simd_double_t quarter = gmx_simd_set1_d(0.25);
2191 const gmx_simd_double_t half = gmx_simd_set1_d(0.5);
2192 const gmx_simd_double_t threequarter = gmx_simd_set1_d(0.75);
2193 gmx_simd_double_t w, q;
2194 gmx_simd_dbool_t m1, m2, m3;
2196 w = gmx_simd_fabs_d(x);
2197 z = gmx_simd_fmadd_d(w, two_over_pi, half);
2198 y = gmx_simd_trunc_d(z);
2199 q = gmx_simd_mul_d(z, quarter);
2200 q = gmx_simd_sub_d(q, gmx_simd_trunc_d(q));
2201 m1 = gmx_simd_cmple_d(quarter, q);
2202 m2 = gmx_simd_cmplt_d(q, half);
2203 m3 = gmx_simd_cmple_d(threequarter, q);
2204 m1 = gmx_simd_and_db(m1, m2);
2205 mask = gmx_simd_or_db(m1, m3);
2206 w = gmx_simd_fnmadd_d(y, argred0, w);
2207 w = gmx_simd_fnmadd_d(y, argred1, w);
2208 w = gmx_simd_fnmadd_d(y, argred2, w);
2209 w = gmx_simd_fnmadd_d(y, argred3, w);
2211 w = gmx_simd_blendv_d(w, gmx_simd_fneg_d(w), mask);
2212 x = gmx_simd_xor_sign_d(w, x);
2214 x2 = gmx_simd_mul_d(x, x);
2215 p = gmx_simd_fmadd_d(CT15, x2, CT14);
2216 p = gmx_simd_fmadd_d(p, x2, CT13);
2217 p = gmx_simd_fmadd_d(p, x2, CT12);
2218 p = gmx_simd_fmadd_d(p, x2, CT11);
2219 p = gmx_simd_fmadd_d(p, x2, CT10);
2220 p = gmx_simd_fmadd_d(p, x2, CT9);
2221 p = gmx_simd_fmadd_d(p, x2, CT8);
2222 p = gmx_simd_fmadd_d(p, x2, CT7);
2223 p = gmx_simd_fmadd_d(p, x2, CT6);
2224 p = gmx_simd_fmadd_d(p, x2, CT5);
2225 p = gmx_simd_fmadd_d(p, x2, CT4);
2226 p = gmx_simd_fmadd_d(p, x2, CT3);
2227 p = gmx_simd_fmadd_d(p, x2, CT2);
2228 p = gmx_simd_fmadd_d(p, x2, CT1);
2229 p = gmx_simd_fmadd_d(x2, gmx_simd_mul_d(p, x), x);
2231 p = gmx_simd_blendv_d( p, gmx_simd_inv_d(p), mask);
2235 /*! \brief SIMD double asin(x).
2237 * \copydetails gmx_simd_asin_f
2239 static gmx_inline gmx_simd_double_t gmx_simdcall
2240 gmx_simd_asin_d(gmx_simd_double_t x)
2242 /* Same algorithm as cephes library */
2243 const gmx_simd_double_t limit1 = gmx_simd_set1_d(0.625);
2244 const gmx_simd_double_t limit2 = gmx_simd_set1_d(1e-8);
2245 const gmx_simd_double_t one = gmx_simd_set1_d(1.0);
2246 const gmx_simd_double_t quarterpi = gmx_simd_set1_d(M_PI/4.0);
2247 const gmx_simd_double_t morebits = gmx_simd_set1_d(6.123233995736765886130e-17);
2249 const gmx_simd_double_t P5 = gmx_simd_set1_d(4.253011369004428248960e-3);
2250 const gmx_simd_double_t P4 = gmx_simd_set1_d(-6.019598008014123785661e-1);
2251 const gmx_simd_double_t P3 = gmx_simd_set1_d(5.444622390564711410273e0);
2252 const gmx_simd_double_t P2 = gmx_simd_set1_d(-1.626247967210700244449e1);
2253 const gmx_simd_double_t P1 = gmx_simd_set1_d(1.956261983317594739197e1);
2254 const gmx_simd_double_t P0 = gmx_simd_set1_d(-8.198089802484824371615e0);
2256 const gmx_simd_double_t Q4 = gmx_simd_set1_d(-1.474091372988853791896e1);
2257 const gmx_simd_double_t Q3 = gmx_simd_set1_d(7.049610280856842141659e1);
2258 const gmx_simd_double_t Q2 = gmx_simd_set1_d(-1.471791292232726029859e2);
2259 const gmx_simd_double_t Q1 = gmx_simd_set1_d(1.395105614657485689735e2);
2260 const gmx_simd_double_t Q0 = gmx_simd_set1_d(-4.918853881490881290097e1);
2262 const gmx_simd_double_t R4 = gmx_simd_set1_d(2.967721961301243206100e-3);
2263 const gmx_simd_double_t R3 = gmx_simd_set1_d(-5.634242780008963776856e-1);
2264 const gmx_simd_double_t R2 = gmx_simd_set1_d(6.968710824104713396794e0);
2265 const gmx_simd_double_t R1 = gmx_simd_set1_d(-2.556901049652824852289e1);
2266 const gmx_simd_double_t R0 = gmx_simd_set1_d(2.853665548261061424989e1);
2268 const gmx_simd_double_t S3 = gmx_simd_set1_d(-2.194779531642920639778e1);
2269 const gmx_simd_double_t S2 = gmx_simd_set1_d(1.470656354026814941758e2);
2270 const gmx_simd_double_t S1 = gmx_simd_set1_d(-3.838770957603691357202e2);
2271 const gmx_simd_double_t S0 = gmx_simd_set1_d(3.424398657913078477438e2);
2273 gmx_simd_double_t xabs;
2274 gmx_simd_double_t zz, ww, z, q, w, zz2, ww2;
2275 gmx_simd_double_t PA, PB;
2276 gmx_simd_double_t QA, QB;
2277 gmx_simd_double_t RA, RB;
2278 gmx_simd_double_t SA, SB;
2279 gmx_simd_double_t nom, denom;
2280 gmx_simd_dbool_t mask;
2282 xabs = gmx_simd_fabs_d(x);
2284 mask = gmx_simd_cmplt_d(limit1, xabs);
2286 zz = gmx_simd_sub_d(one, xabs);
2287 ww = gmx_simd_mul_d(xabs, xabs);
2288 zz2 = gmx_simd_mul_d(zz, zz);
2289 ww2 = gmx_simd_mul_d(ww, ww);
2292 RA = gmx_simd_mul_d(R4, zz2);
2293 RB = gmx_simd_mul_d(R3, zz2);
2294 RA = gmx_simd_add_d(RA, R2);
2295 RB = gmx_simd_add_d(RB, R1);
2296 RA = gmx_simd_mul_d(RA, zz2);
2297 RB = gmx_simd_mul_d(RB, zz);
2298 RA = gmx_simd_add_d(RA, R0);
2299 RA = gmx_simd_add_d(RA, RB);
2302 SB = gmx_simd_mul_d(S3, zz2);
2303 SA = gmx_simd_add_d(zz2, S2);
2304 SB = gmx_simd_add_d(SB, S1);
2305 SA = gmx_simd_mul_d(SA, zz2);
2306 SB = gmx_simd_mul_d(SB, zz);
2307 SA = gmx_simd_add_d(SA, S0);
2308 SA = gmx_simd_add_d(SA, SB);
2311 PA = gmx_simd_mul_d(P5, ww2);
2312 PB = gmx_simd_mul_d(P4, ww2);
2313 PA = gmx_simd_add_d(PA, P3);
2314 PB = gmx_simd_add_d(PB, P2);
2315 PA = gmx_simd_mul_d(PA, ww2);
2316 PB = gmx_simd_mul_d(PB, ww2);
2317 PA = gmx_simd_add_d(PA, P1);
2318 PB = gmx_simd_add_d(PB, P0);
2319 PA = gmx_simd_mul_d(PA, ww);
2320 PA = gmx_simd_add_d(PA, PB);
2323 QB = gmx_simd_mul_d(Q4, ww2);
2324 QA = gmx_simd_add_d(ww2, Q3);
2325 QB = gmx_simd_add_d(QB, Q2);
2326 QA = gmx_simd_mul_d(QA, ww2);
2327 QB = gmx_simd_mul_d(QB, ww2);
2328 QA = gmx_simd_add_d(QA, Q1);
2329 QB = gmx_simd_add_d(QB, Q0);
2330 QA = gmx_simd_mul_d(QA, ww);
2331 QA = gmx_simd_add_d(QA, QB);
2333 RA = gmx_simd_mul_d(RA, zz);
2334 PA = gmx_simd_mul_d(PA, ww);
2336 nom = gmx_simd_blendv_d( PA, RA, mask );
2337 denom = gmx_simd_blendv_d( QA, SA, mask );
2339 q = gmx_simd_mul_d( nom, gmx_simd_inv_d(denom) );
2341 zz = gmx_simd_add_d(zz, zz);
2342 zz = gmx_simd_sqrt_d(zz);
2343 z = gmx_simd_sub_d(quarterpi, zz);
2344 zz = gmx_simd_mul_d(zz, q);
2345 zz = gmx_simd_sub_d(zz, morebits);
2346 z = gmx_simd_sub_d(z, zz);
2347 z = gmx_simd_add_d(z, quarterpi);
2349 w = gmx_simd_mul_d(xabs, q);
2350 w = gmx_simd_add_d(w, xabs);
2352 z = gmx_simd_blendv_d( w, z, mask );
2354 mask = gmx_simd_cmplt_d(limit2, xabs);
2355 z = gmx_simd_blendv_d( xabs, z, mask );
2357 z = gmx_simd_xor_sign_d(z, x);
2362 /*! \brief SIMD double acos(x).
2364 * \copydetails gmx_simd_acos_f
2366 static gmx_inline gmx_simd_double_t gmx_simdcall
2367 gmx_simd_acos_d(gmx_simd_double_t x)
2369 const gmx_simd_double_t one = gmx_simd_set1_d(1.0);
2370 const gmx_simd_double_t half = gmx_simd_set1_d(0.5);
2371 const gmx_simd_double_t quarterpi0 = gmx_simd_set1_d(7.85398163397448309616e-1);
2372 const gmx_simd_double_t quarterpi1 = gmx_simd_set1_d(6.123233995736765886130e-17);
2374 gmx_simd_dbool_t mask1;
2375 gmx_simd_double_t z, z1, z2;
2377 mask1 = gmx_simd_cmplt_d(half, x);
2378 z1 = gmx_simd_mul_d(half, gmx_simd_sub_d(one, x));
2379 z1 = gmx_simd_sqrt_d(z1);
2380 z = gmx_simd_blendv_d( x, z1, mask1 );
2382 z = gmx_simd_asin_d(z);
2384 z1 = gmx_simd_add_d(z, z);
2386 z2 = gmx_simd_sub_d(quarterpi0, z);
2387 z2 = gmx_simd_add_d(z2, quarterpi1);
2388 z2 = gmx_simd_add_d(z2, quarterpi0);
2390 z = gmx_simd_blendv_d(z2, z1, mask1);
2395 /*! \brief SIMD double atan(x).
2397 * \copydetails gmx_simd_atan_f
2399 static gmx_inline gmx_simd_double_t gmx_simdcall
2400 gmx_simd_atan_d(gmx_simd_double_t x)
2402 /* Same algorithm as cephes library */
2403 const gmx_simd_double_t limit1 = gmx_simd_set1_d(0.66);
2404 const gmx_simd_double_t limit2 = gmx_simd_set1_d(2.41421356237309504880);
2405 const gmx_simd_double_t quarterpi = gmx_simd_set1_d(M_PI/4.0);
2406 const gmx_simd_double_t halfpi = gmx_simd_set1_d(M_PI/2.0);
2407 const gmx_simd_double_t mone = gmx_simd_set1_d(-1.0);
2408 const gmx_simd_double_t morebits1 = gmx_simd_set1_d(0.5*6.123233995736765886130E-17);
2409 const gmx_simd_double_t morebits2 = gmx_simd_set1_d(6.123233995736765886130E-17);
2411 const gmx_simd_double_t P4 = gmx_simd_set1_d(-8.750608600031904122785E-1);
2412 const gmx_simd_double_t P3 = gmx_simd_set1_d(-1.615753718733365076637E1);
2413 const gmx_simd_double_t P2 = gmx_simd_set1_d(-7.500855792314704667340E1);
2414 const gmx_simd_double_t P1 = gmx_simd_set1_d(-1.228866684490136173410E2);
2415 const gmx_simd_double_t P0 = gmx_simd_set1_d(-6.485021904942025371773E1);
2417 const gmx_simd_double_t Q4 = gmx_simd_set1_d(2.485846490142306297962E1);
2418 const gmx_simd_double_t Q3 = gmx_simd_set1_d(1.650270098316988542046E2);
2419 const gmx_simd_double_t Q2 = gmx_simd_set1_d(4.328810604912902668951E2);
2420 const gmx_simd_double_t Q1 = gmx_simd_set1_d(4.853903996359136964868E2);
2421 const gmx_simd_double_t Q0 = gmx_simd_set1_d(1.945506571482613964425E2);
2423 gmx_simd_double_t y, xabs, t1, t2;
2424 gmx_simd_double_t z, z2;
2425 gmx_simd_double_t P_A, P_B, Q_A, Q_B;
2426 gmx_simd_dbool_t mask1, mask2;
2428 xabs = gmx_simd_fabs_d(x);
2430 mask1 = gmx_simd_cmplt_d(limit1, xabs);
2431 mask2 = gmx_simd_cmplt_d(limit2, xabs);
2433 t1 = gmx_simd_mul_d(gmx_simd_add_d(xabs, mone), gmx_simd_inv_d(gmx_simd_sub_d(xabs, mone)));
2434 t2 = gmx_simd_mul_d(mone, gmx_simd_inv_d(xabs));
2436 y = gmx_simd_blendzero_d(quarterpi, mask1);
2437 y = gmx_simd_blendv_d(y, halfpi, mask2);
2438 xabs = gmx_simd_blendv_d(xabs, t1, mask1);
2439 xabs = gmx_simd_blendv_d(xabs, t2, mask2);
2441 z = gmx_simd_mul_d(xabs, xabs);
2442 z2 = gmx_simd_mul_d(z, z);
2444 P_A = gmx_simd_mul_d(P4, z2);
2445 P_B = gmx_simd_mul_d(P3, z2);
2446 P_A = gmx_simd_add_d(P_A, P2);
2447 P_B = gmx_simd_add_d(P_B, P1);
2448 P_A = gmx_simd_mul_d(P_A, z2);
2449 P_B = gmx_simd_mul_d(P_B, z);
2450 P_A = gmx_simd_add_d(P_A, P0);
2451 P_A = gmx_simd_add_d(P_A, P_B);
2454 Q_B = gmx_simd_mul_d(Q4, z2);
2455 Q_A = gmx_simd_add_d(z2, Q3);
2456 Q_B = gmx_simd_add_d(Q_B, Q2);
2457 Q_A = gmx_simd_mul_d(Q_A, z2);
2458 Q_B = gmx_simd_mul_d(Q_B, z2);
2459 Q_A = gmx_simd_add_d(Q_A, Q1);
2460 Q_B = gmx_simd_add_d(Q_B, Q0);
2461 Q_A = gmx_simd_mul_d(Q_A, z);
2462 Q_A = gmx_simd_add_d(Q_A, Q_B);
2464 z = gmx_simd_mul_d(z, P_A);
2465 z = gmx_simd_mul_d(z, gmx_simd_inv_d(Q_A));
2466 z = gmx_simd_mul_d(z, xabs);
2467 z = gmx_simd_add_d(z, xabs);
2469 t1 = gmx_simd_blendzero_d(morebits1, mask1);
2470 t1 = gmx_simd_blendv_d(t1, morebits2, mask2);
2472 z = gmx_simd_add_d(z, t1);
2473 y = gmx_simd_add_d(y, z);
2475 y = gmx_simd_xor_sign_d(y, x);
2480 /*! \brief SIMD double atan2(y,x).
2482 * \copydetails gmx_simd_atan2_f
2484 static gmx_inline gmx_simd_double_t gmx_simdcall
2485 gmx_simd_atan2_d(gmx_simd_double_t y, gmx_simd_double_t x)
2487 const gmx_simd_double_t pi = gmx_simd_set1_d(M_PI);
2488 const gmx_simd_double_t halfpi = gmx_simd_set1_d(M_PI/2.0);
2489 gmx_simd_double_t xinv, p, aoffset;
2490 gmx_simd_dbool_t mask_x0, mask_y0, mask_xlt0, mask_ylt0;
2492 mask_x0 = gmx_simd_cmpeq_d(x, gmx_simd_setzero_d());
2493 mask_y0 = gmx_simd_cmpeq_d(y, gmx_simd_setzero_d());
2494 mask_xlt0 = gmx_simd_cmplt_d(x, gmx_simd_setzero_d());
2495 mask_ylt0 = gmx_simd_cmplt_d(y, gmx_simd_setzero_d());
2497 aoffset = gmx_simd_blendzero_d(halfpi, mask_x0);
2498 aoffset = gmx_simd_blendnotzero_d(aoffset, mask_y0);
2500 aoffset = gmx_simd_blendv_d(aoffset, pi, mask_xlt0);
2501 aoffset = gmx_simd_blendv_d(aoffset, gmx_simd_fneg_d(aoffset), mask_ylt0);
2503 xinv = gmx_simd_blendnotzero_d(gmx_simd_inv_d(x), mask_x0);
2504 p = gmx_simd_mul_d(y, xinv);
2505 p = gmx_simd_atan_d(p);
2506 p = gmx_simd_add_d(p, aoffset);
2512 /*! \brief Calculate the force correction due to PME analytically for SIMD double.
2514 * \copydetails gmx_simd_pmecorrF_f
2516 static gmx_inline gmx_simd_double_t gmx_simdcall
2517 gmx_simd_pmecorrF_d(gmx_simd_double_t z2)
2519 const gmx_simd_double_t FN10 = gmx_simd_set1_d(-8.0072854618360083154e-14);
2520 const gmx_simd_double_t FN9 = gmx_simd_set1_d(1.1859116242260148027e-11);
2521 const gmx_simd_double_t FN8 = gmx_simd_set1_d(-8.1490406329798423616e-10);
2522 const gmx_simd_double_t FN7 = gmx_simd_set1_d(3.4404793543907847655e-8);
2523 const gmx_simd_double_t FN6 = gmx_simd_set1_d(-9.9471420832602741006e-7);
2524 const gmx_simd_double_t FN5 = gmx_simd_set1_d(0.000020740315999115847456);
2525 const gmx_simd_double_t FN4 = gmx_simd_set1_d(-0.00031991745139313364005);
2526 const gmx_simd_double_t FN3 = gmx_simd_set1_d(0.0035074449373659008203);
2527 const gmx_simd_double_t FN2 = gmx_simd_set1_d(-0.031750380176100813405);
2528 const gmx_simd_double_t FN1 = gmx_simd_set1_d(0.13884101728898463426);
2529 const gmx_simd_double_t FN0 = gmx_simd_set1_d(-0.75225277815249618847);
2531 const gmx_simd_double_t FD5 = gmx_simd_set1_d(0.000016009278224355026701);
2532 const gmx_simd_double_t FD4 = gmx_simd_set1_d(0.00051055686934806966046);
2533 const gmx_simd_double_t FD3 = gmx_simd_set1_d(0.0081803507497974289008);
2534 const gmx_simd_double_t FD2 = gmx_simd_set1_d(0.077181146026670287235);
2535 const gmx_simd_double_t FD1 = gmx_simd_set1_d(0.41543303143712535988);
2536 const gmx_simd_double_t FD0 = gmx_simd_set1_d(1.0);
2538 gmx_simd_double_t z4;
2539 gmx_simd_double_t polyFN0, polyFN1, polyFD0, polyFD1;
2541 z4 = gmx_simd_mul_d(z2, z2);
2543 polyFD1 = gmx_simd_fmadd_d(FD5, z4, FD3);
2544 polyFD1 = gmx_simd_fmadd_d(polyFD1, z4, FD1);
2545 polyFD1 = gmx_simd_mul_d(polyFD1, z2);
2546 polyFD0 = gmx_simd_fmadd_d(FD4, z4, FD2);
2547 polyFD0 = gmx_simd_fmadd_d(polyFD0, z4, FD0);
2548 polyFD0 = gmx_simd_add_d(polyFD0, polyFD1);
2550 polyFD0 = gmx_simd_inv_d(polyFD0);
2552 polyFN0 = gmx_simd_fmadd_d(FN10, z4, FN8);
2553 polyFN0 = gmx_simd_fmadd_d(polyFN0, z4, FN6);
2554 polyFN0 = gmx_simd_fmadd_d(polyFN0, z4, FN4);
2555 polyFN0 = gmx_simd_fmadd_d(polyFN0, z4, FN2);
2556 polyFN0 = gmx_simd_fmadd_d(polyFN0, z4, FN0);
2557 polyFN1 = gmx_simd_fmadd_d(FN9, z4, FN7);
2558 polyFN1 = gmx_simd_fmadd_d(polyFN1, z4, FN5);
2559 polyFN1 = gmx_simd_fmadd_d(polyFN1, z4, FN3);
2560 polyFN1 = gmx_simd_fmadd_d(polyFN1, z4, FN1);
2561 polyFN0 = gmx_simd_fmadd_d(polyFN1, z2, polyFN0);
2564 return gmx_simd_mul_d(polyFN0, polyFD0);
2569 /*! \brief Calculate the potential correction due to PME analytically for SIMD double.
2571 * \copydetails gmx_simd_pmecorrV_f
2573 static gmx_inline gmx_simd_double_t gmx_simdcall
2574 gmx_simd_pmecorrV_d(gmx_simd_double_t z2)
2576 const gmx_simd_double_t VN9 = gmx_simd_set1_d(-9.3723776169321855475e-13);
2577 const gmx_simd_double_t VN8 = gmx_simd_set1_d(1.2280156762674215741e-10);
2578 const gmx_simd_double_t VN7 = gmx_simd_set1_d(-7.3562157912251309487e-9);
2579 const gmx_simd_double_t VN6 = gmx_simd_set1_d(2.6215886208032517509e-7);
2580 const gmx_simd_double_t VN5 = gmx_simd_set1_d(-4.9532491651265819499e-6);
2581 const gmx_simd_double_t VN4 = gmx_simd_set1_d(0.00025907400778966060389);
2582 const gmx_simd_double_t VN3 = gmx_simd_set1_d(0.0010585044856156469792);
2583 const gmx_simd_double_t VN2 = gmx_simd_set1_d(0.045247661136833092885);
2584 const gmx_simd_double_t VN1 = gmx_simd_set1_d(0.11643931522926034421);
2585 const gmx_simd_double_t VN0 = gmx_simd_set1_d(1.1283791671726767970);
2587 const gmx_simd_double_t VD5 = gmx_simd_set1_d(0.000021784709867336150342);
2588 const gmx_simd_double_t VD4 = gmx_simd_set1_d(0.00064293662010911388448);
2589 const gmx_simd_double_t VD3 = gmx_simd_set1_d(0.0096311444822588683504);
2590 const gmx_simd_double_t VD2 = gmx_simd_set1_d(0.085608012351550627051);
2591 const gmx_simd_double_t VD1 = gmx_simd_set1_d(0.43652499166614811084);
2592 const gmx_simd_double_t VD0 = gmx_simd_set1_d(1.0);
2594 gmx_simd_double_t z4;
2595 gmx_simd_double_t polyVN0, polyVN1, polyVD0, polyVD1;
2597 z4 = gmx_simd_mul_d(z2, z2);
2599 polyVD1 = gmx_simd_fmadd_d(VD5, z4, VD3);
2600 polyVD0 = gmx_simd_fmadd_d(VD4, z4, VD2);
2601 polyVD1 = gmx_simd_fmadd_d(polyVD1, z4, VD1);
2602 polyVD0 = gmx_simd_fmadd_d(polyVD0, z4, VD0);
2603 polyVD0 = gmx_simd_fmadd_d(polyVD1, z2, polyVD0);
2605 polyVD0 = gmx_simd_inv_d(polyVD0);
2607 polyVN1 = gmx_simd_fmadd_d(VN9, z4, VN7);
2608 polyVN0 = gmx_simd_fmadd_d(VN8, z4, VN6);
2609 polyVN1 = gmx_simd_fmadd_d(polyVN1, z4, VN5);
2610 polyVN0 = gmx_simd_fmadd_d(polyVN0, z4, VN4);
2611 polyVN1 = gmx_simd_fmadd_d(polyVN1, z4, VN3);
2612 polyVN0 = gmx_simd_fmadd_d(polyVN0, z4, VN2);
2613 polyVN1 = gmx_simd_fmadd_d(polyVN1, z4, VN1);
2614 polyVN0 = gmx_simd_fmadd_d(polyVN0, z4, VN0);
2615 polyVN0 = gmx_simd_fmadd_d(polyVN1, z2, polyVN0);
2617 return gmx_simd_mul_d(polyVN0, polyVD0);
2625 /*! \name SIMD4 math functions
2627 * \note Only a subset of the math functions are implemented for SIMD4.
2632 #ifdef GMX_SIMD4_HAVE_FLOAT
2634 /*************************************************************************
2635 * SINGLE PRECISION SIMD4 MATH FUNCTIONS - JUST A SMALL SUBSET SUPPORTED *
2636 *************************************************************************/
2638 /*! \brief SIMD4 utility function to sum a+b+c+d for SIMD4 floats.
2640 * \copydetails gmx_simd_sum4_f
2642 static gmx_inline gmx_simd4_float_t gmx_simdcall
2643 gmx_simd4_sum4_f(gmx_simd4_float_t a, gmx_simd4_float_t b,
2644 gmx_simd4_float_t c, gmx_simd4_float_t d)
2646 return gmx_simd4_add_f(gmx_simd4_add_f(a, b), gmx_simd4_add_f(c, d));
2649 /*! \brief Perform one Newton-Raphson iteration to improve 1/sqrt(x) for SIMD4 float.
2651 * \copydetails gmx_simd_rsqrt_iter_f
2653 static gmx_inline gmx_simd4_float_t gmx_simdcall
2654 gmx_simd4_rsqrt_iter_f(gmx_simd4_float_t lu, gmx_simd4_float_t x)
2656 # ifdef GMX_SIMD_HAVE_FMA
2657 return gmx_simd4_fmadd_f(gmx_simd4_fnmadd_f(x, gmx_simd4_mul_f(lu, lu), gmx_simd4_set1_f(1.0f)), gmx_simd4_mul_f(lu, gmx_simd4_set1_f(0.5f)), lu);
2659 return gmx_simd4_mul_f(gmx_simd4_set1_f(0.5f), gmx_simd4_mul_f(gmx_simd4_sub_f(gmx_simd4_set1_f(3.0f), gmx_simd4_mul_f(gmx_simd4_mul_f(lu, lu), x)), lu));
2663 /*! \brief Calculate 1/sqrt(x) for SIMD4 float.
2665 * \copydetails gmx_simd_invsqrt_f
2667 static gmx_inline gmx_simd4_float_t gmx_simdcall
2668 gmx_simd4_invsqrt_f(gmx_simd4_float_t x)
2670 gmx_simd4_float_t lu = gmx_simd4_rsqrt_f(x);
2671 #if (GMX_SIMD_RSQRT_BITS < GMX_SIMD_MATH_TARGET_SINGLE_BITS)
2672 lu = gmx_simd4_rsqrt_iter_f(lu, x);
2674 #if (GMX_SIMD_RSQRT_BITS*2 < GMX_SIMD_MATH_TARGET_SINGLE_BITS)
2675 lu = gmx_simd4_rsqrt_iter_f(lu, x);
2677 #if (GMX_SIMD_RSQRT_BITS*4 < GMX_SIMD_MATH_TARGET_SINGLE_BITS)
2678 lu = gmx_simd4_rsqrt_iter_f(lu, x);
2683 #endif /* GMX_SIMD4_HAVE_FLOAT */
2687 #ifdef GMX_SIMD4_HAVE_DOUBLE
2688 /*************************************************************************
2689 * DOUBLE PRECISION SIMD4 MATH FUNCTIONS - JUST A SMALL SUBSET SUPPORTED *
2690 *************************************************************************/
2693 /*! \brief SIMD4 utility function to sum a+b+c+d for SIMD4 doubles.
2695 * \copydetails gmx_simd_sum4_f
2697 static gmx_inline gmx_simd4_double_t gmx_simdcall
2698 gmx_simd4_sum4_d(gmx_simd4_double_t a, gmx_simd4_double_t b,
2699 gmx_simd4_double_t c, gmx_simd4_double_t d)
2701 return gmx_simd4_add_d(gmx_simd4_add_d(a, b), gmx_simd4_add_d(c, d));
2704 /*! \brief Perform one Newton-Raphson iteration to improve 1/sqrt(x) for SIMD4 double.
2706 * \copydetails gmx_simd_rsqrt_iter_f
2708 static gmx_inline gmx_simd4_double_t gmx_simdcall
2709 gmx_simd4_rsqrt_iter_d(gmx_simd4_double_t lu, gmx_simd4_double_t x)
2711 #ifdef GMX_SIMD_HAVE_FMA
2712 return gmx_simd4_fmadd_d(gmx_simd4_fnmadd_d(x, gmx_simd4_mul_d(lu, lu), gmx_simd4_set1_d(1.0)), gmx_simd4_mul_d(lu, gmx_simd4_set1_d(0.5)), lu);
2714 return gmx_simd4_mul_d(gmx_simd4_set1_d(0.5), gmx_simd4_mul_d(gmx_simd4_sub_d(gmx_simd4_set1_d(3.0), gmx_simd4_mul_d(gmx_simd4_mul_d(lu, lu), x)), lu));
2718 /*! \brief Calculate 1/sqrt(x) for SIMD4 double.
2720 * \copydetails gmx_simd_invsqrt_f
2722 static gmx_inline gmx_simd4_double_t gmx_simdcall
2723 gmx_simd4_invsqrt_d(gmx_simd4_double_t x)
2725 gmx_simd4_double_t lu = gmx_simd4_rsqrt_d(x);
2726 #if (GMX_SIMD_RSQRT_BITS < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
2727 lu = gmx_simd4_rsqrt_iter_d(lu, x);
2729 #if (GMX_SIMD_RSQRT_BITS*2 < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
2730 lu = gmx_simd4_rsqrt_iter_d(lu, x);
2732 #if (GMX_SIMD_RSQRT_BITS*4 < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
2733 lu = gmx_simd4_rsqrt_iter_d(lu, x);
2735 #if (GMX_SIMD_RSQRT_BITS*8 < GMX_SIMD_MATH_TARGET_DOUBLE_BITS)
2736 lu = gmx_simd4_rsqrt_iter_d(lu, x);
2740 #endif /* GMX_SIMD4_HAVE_DOUBLE */
2745 /* Set defines based on default Gromacs precision */
2747 /* Documentation in single branch below */
2748 # define gmx_simd_sum4_r gmx_simd_sum4_d
2749 # define gmx_simd_xor_sign_r gmx_simd_xor_sign_d
2750 # define gmx_simd_invsqrt_r gmx_simd_invsqrt_d
2751 # define gmx_simd_invsqrt_pair_r gmx_simd_invsqrt_pair_d
2752 # define gmx_simd_sqrt_r gmx_simd_sqrt_d
2753 # define gmx_simd_inv_r gmx_simd_inv_d
2754 # define gmx_simd_log_r gmx_simd_log_d
2755 # define gmx_simd_exp2_r gmx_simd_exp2_d
2756 # define gmx_simd_exp_r gmx_simd_exp_d
2757 # define gmx_simd_erf_r gmx_simd_erf_d
2758 # define gmx_simd_erfc_r gmx_simd_erfc_d
2759 # define gmx_simd_sincos_r gmx_simd_sincos_d
2760 # define gmx_simd_sin_r gmx_simd_sin_d
2761 # define gmx_simd_cos_r gmx_simd_cos_d
2762 # define gmx_simd_tan_r gmx_simd_tan_d
2763 # define gmx_simd_asin_r gmx_simd_asin_d
2764 # define gmx_simd_acos_r gmx_simd_acos_d
2765 # define gmx_simd_atan_r gmx_simd_atan_d
2766 # define gmx_simd_atan2_r gmx_simd_atan2_d
2767 # define gmx_simd_pmecorrF_r gmx_simd_pmecorrF_d
2768 # define gmx_simd_pmecorrV_r gmx_simd_pmecorrV_d
2769 # define gmx_simd4_sum4_r gmx_simd4_sum4_d
2770 # define gmx_simd4_invsqrt_r gmx_simd4_invsqrt_d
2772 #else /* GMX_DOUBLE */
2774 /*! \name Real-precision SIMD math functions
2776 * These are the ones you should typically call in Gromacs.
2780 /*! \brief SIMD utility function to sum a+b+c+d for SIMD reals.
2782 * \copydetails gmx_simd_sum4_f
2784 # define gmx_simd_sum4_r gmx_simd_sum4_f
2786 /*! \brief Return -a if b is negative, SIMD real.
2788 * \copydetails gmx_simd_xor_sign_f
2790 # define gmx_simd_xor_sign_r gmx_simd_xor_sign_f
2792 /*! \brief Calculate 1/sqrt(x) for SIMD real.
2794 * \copydetails gmx_simd_invsqrt_f
2796 # define gmx_simd_invsqrt_r gmx_simd_invsqrt_f
2798 /*! \brief Calculate 1/sqrt(x) for two SIMD reals.
2800 * \copydetails gmx_simd_invsqrt_pair_f
2802 # define gmx_simd_invsqrt_pair_r gmx_simd_invsqrt_pair_f
2804 /*! \brief Calculate sqrt(x) correctly for SIMD real, including argument 0.0.
2806 * \copydetails gmx_simd_sqrt_f
2808 # define gmx_simd_sqrt_r gmx_simd_sqrt_f
2810 /*! \brief Calculate 1/x for SIMD real.
2812 * \copydetails gmx_simd_inv_f
2814 # define gmx_simd_inv_r gmx_simd_inv_f
2816 /*! \brief SIMD real log(x). This is the natural logarithm.
2818 * \copydetails gmx_simd_log_f
2820 # define gmx_simd_log_r gmx_simd_log_f
2822 /*! \brief SIMD real 2^x.
2824 * \copydetails gmx_simd_exp2_f
2826 # define gmx_simd_exp2_r gmx_simd_exp2_f
2828 /*! \brief SIMD real e^x.
2830 * \copydetails gmx_simd_exp_f
2832 # define gmx_simd_exp_r gmx_simd_exp_f
2834 /*! \brief SIMD real erf(x).
2836 * \copydetails gmx_simd_erf_f
2838 # define gmx_simd_erf_r gmx_simd_erf_f
2840 /*! \brief SIMD real erfc(x).
2842 * \copydetails gmx_simd_erfc_f
2844 # define gmx_simd_erfc_r gmx_simd_erfc_f
2846 /*! \brief SIMD real sin \& cos.
2848 * \copydetails gmx_simd_sincos_f
2850 # define gmx_simd_sincos_r gmx_simd_sincos_f
2852 /*! \brief SIMD real sin(x).
2854 * \copydetails gmx_simd_sin_f
2856 # define gmx_simd_sin_r gmx_simd_sin_f
2858 /*! \brief SIMD real cos(x).
2860 * \copydetails gmx_simd_cos_f
2862 # define gmx_simd_cos_r gmx_simd_cos_f
2864 /*! \brief SIMD real tan(x).
2866 * \copydetails gmx_simd_tan_f
2868 # define gmx_simd_tan_r gmx_simd_tan_f
2870 /*! \brief SIMD real asin(x).
2872 * \copydetails gmx_simd_asin_f
2874 # define gmx_simd_asin_r gmx_simd_asin_f
2876 /*! \brief SIMD real acos(x).
2878 * \copydetails gmx_simd_acos_f
2880 # define gmx_simd_acos_r gmx_simd_acos_f
2882 /*! \brief SIMD real atan(x).
2884 * \copydetails gmx_simd_atan_f
2886 # define gmx_simd_atan_r gmx_simd_atan_f
2888 /*! \brief SIMD real atan2(y,x).
2890 * \copydetails gmx_simd_atan2_f
2892 # define gmx_simd_atan2_r gmx_simd_atan2_f
2894 /*! \brief SIMD Analytic PME force correction.
2896 * \copydetails gmx_simd_pmecorrF_f
2898 # define gmx_simd_pmecorrF_r gmx_simd_pmecorrF_f
2900 /*! \brief SIMD Analytic PME potential correction.
2902 * \copydetails gmx_simd_pmecorrV_f
2904 # define gmx_simd_pmecorrV_r gmx_simd_pmecorrV_f
2907 * \name SIMD4 math functions
2911 /*! \brief SIMD4 utility function to sum a+b+c+d for SIMD4 reals.
2913 * \copydetails gmx_simd_sum4_f
2915 # define gmx_simd4_sum4_r gmx_simd4_sum4_f
2917 /*! \brief Calculate 1/sqrt(x) for SIMD4 real.
2919 * \copydetails gmx_simd_invsqrt_f
2921 # define gmx_simd4_invsqrt_r gmx_simd4_invsqrt_f
2925 #endif /* GMX_DOUBLE */
2930 #endif /* GMX_SIMD_SIMD_MATH_H_ */