Add Power8 VSX SIMD support
[alexxy/gromacs.git] / src / gromacs / simd / simd.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 /*! \libinternal
37  * \defgroup module_simd SIMD intrinsics interface (simd)
38  * \ingroup group_utilitymodules
39  *
40  * \brief Provides an architecture-independent way of doing SIMD coding.
41  *
42  * Overview of the SIMD implementation is provided in \ref page_simd.
43  * The details are documented in simd.h and the reference implementation
44  * impl_reference.h.
45  *
46  * \author Erik Lindahl <erik.lindahl@scilifelab.se>
47  */
48
49 #ifndef GMX_SIMD_SIMD_H
50 #define GMX_SIMD_SIMD_H
51
52 /*! \libinternal \file
53  *
54  * \brief Definitions, capabilities, and wrappers for SIMD module.
55  *
56  * The macros in this file are intended to be used for writing
57  * architecture-independent SIMD intrinsics code.
58  * To support a new architecture, adding a new sub-include with macros here
59  * should be (nearly) all that is needed.
60  *
61  * The defines in this top-level file will set default Gromacs real precision
62  * operations to either single or double precision based on whether
63  * GMX_DOUBLE is defined. The actual implementation - including e.g.
64  * conversion operations specifically between single and double - is documented
65  * in impl_reference.h.
66  *
67  * \author Erik Lindahl <erik.lindahl@scilifelab.se>
68  *
69  * \inlibraryapi
70  * \ingroup module_simd
71  */
72
73 #include "config.h"
74
75 #include <stddef.h>
76
77 #include "gromacs/utility/basedefinitions.h"
78
79 /* Forward declarations so memory allocation can be used in implementations */
80 static gmx_inline float *  gmx_simd_align_f(float *p);
81 static gmx_inline double * gmx_simd_align_d(double *p);
82 static gmx_inline int *    gmx_simd_align_fi(int *p);
83 static gmx_inline int *    gmx_simd_align_di(int *p);
84 static gmx_inline float *  gmx_simd4_align_f(float *p);
85 static gmx_inline double * gmx_simd4_align_d(double *p);
86
87 /*! \cond libapi */
88 /*! \addtogroup module_simd */
89 /*! \{ */
90
91 /*! \name SIMD predefined macros to describe high-level capabilities
92  *
93  *  These macros are used to describe the features available in default
94  *  Gromacs real precision. They are set from the lower-level implementation
95  *  files that have macros describing single and double precision individually,
96  *  as well as the implementation details.
97  *  \{
98  */
99
100 /*! \brief
101  *  GMX_SIMD indicates that some sort of SIMD support is present in software.
102  *
103  * It is disabled if no architecture, neither reference SIMD, has been selected.
104  */
105 #define GMX_SIMD
106
107
108 /* Intel MIC is a bit special since it is a co-processor. This means the rest
109  * of GROMACS (which runs on the CPU) should use a default SIMD set like AVX,
110  * while the part running on the coprocessor defines __MIC__. All functions in
111  * this SIMD module are static, so it will work perfectly fine to include this
112  * file with different SIMD definitions for different files.
113  */
114 #if defined __MIC__
115 #    include "impl_intel_mic/impl_intel_mic.h"
116 #elif defined GMX_SIMD_X86_AVX2_256
117 #    include "impl_x86_avx2_256/impl_x86_avx2_256.h"
118 #elif defined GMX_SIMD_X86_AVX_256
119 #    include "impl_x86_avx_256/impl_x86_avx_256.h"
120 #elif defined GMX_SIMD_X86_AVX_128_FMA
121 #    include "impl_x86_avx_128_fma/impl_x86_avx_128_fma.h"
122 #elif defined GMX_SIMD_X86_SSE4_1
123 #    include "impl_x86_sse4_1/impl_x86_sse4_1.h"
124 #elif defined GMX_SIMD_X86_SSE2
125 #    include "impl_x86_sse2/impl_x86_sse2.h"
126 #elif defined GMX_SIMD_ARM_NEON
127 #    include "impl_arm_neon/impl_arm_neon.h"
128 #elif defined GMX_SIMD_ARM_NEON_ASIMD
129 #    include "impl_arm_neon_asimd/impl_arm_neon_asimd.h"
130 #elif defined GMX_SIMD_IBM_QPX
131 #    include "impl_ibm_qpx/impl_ibm_qpx.h"
132 #elif defined GMX_SIMD_IBM_VMX
133 #    include "impl_ibm_vmx/impl_ibm_vmx.h"
134 #elif defined GMX_SIMD_IBM_VSX
135 #    include "impl_ibm_vsx/impl_ibm_vsx.h"
136 #elif defined GMX_SIMD_SPARC64_HPC_ACE
137 #    include "impl_sparc64_hpc_ace/impl_sparc64_hpc_ace.h"
138 #elif (defined GMX_SIMD_REFERENCE) || (defined DOXYGEN)
139 /* Plain C SIMD reference implementation, also serves as documentation.
140  * For now this code path will also be taken for Sparc64_HPC_ACE since we have
141  * not yet added the verlet kernel extensions there. The group kernels do not
142  * depend on this file, so they will still be accelerated with SIMD.
143  */
144 #    include "impl_reference/impl_reference.h"
145 #else
146 /* Turn off the GMX_SIMD flag if we do not even have reference support */
147 #    undef GMX_SIMD
148 #endif
149
150 /*! \brief
151  * SIMD4 width is always 4, but use this for clarity in definitions.
152  *
153  * It improves code readability to allocate e.g. 2*GMX_SIMD4_WIDTH instead of 8.
154  */
155 #define GMX_SIMD4_WIDTH    4
156
157 /*! \} */
158
159 /*! \name SIMD memory alignment operations
160  *  \{
161  */
162
163 /*! \brief
164  * Align a float pointer for usage with SIMD instructions.
165  *
166  * You should typically \a not call this function directly (unless you explicitly
167  * want single precision even when GMX_DOUBLE is set), but use the
168  * \ref gmx_simd_align_r macro to align memory in default Gromacs real precision.
169  *
170  * \param  p Pointer to memory, allocate at least \ref GMX_SIMD_FLOAT_WIDTH extra elements.
171  *
172  * \return Aligned pointer (>=p) suitable for loading/storing float fp SIMD.
173  *         If \ref GMX_SIMD_HAVE_FLOAT is not set, p will be returned unchanged.
174  *
175  * Start by allocating an extra \ref GMX_SIMD_FLOAT_WIDTH float elements of memory,
176  * and then call this function. The returned pointer will be greater or equal
177  * to the one you provided, and point to an address inside your provided memory
178  * that is aligned to the SIMD width.
179  */
180 static gmx_inline float *
181 gmx_simd_align_f(float *p)
182 {
183 #    ifdef GMX_SIMD_HAVE_FLOAT
184     return (float *)(((size_t)((p)+GMX_SIMD_FLOAT_WIDTH-1)) & (~((size_t)(GMX_SIMD_FLOAT_WIDTH*sizeof(float)-1))));
185 #    else
186     return p;
187 #    endif
188 }
189
190 /*!  \brief
191  * Align a double pointer for usage with SIMD instructions.
192  *
193  * You should typically \a not call this function directly (unless you explicitly
194  * want double precision even when GMX_DOUBLE is not set), but use the
195  * \ref gmx_simd_align_r macro to align memory in default Gromacs real precision.
196  *
197  * \param  p Pointer to memory, allocate at least \ref GMX_SIMD_DOUBLE_WIDTH extra elements.
198  *
199  * \return Aligned pointer (>=p) suitable for loading/storing double fp SIMD.
200  *         If \ref GMX_SIMD_HAVE_DOUBLE is not set, p will be returned unchanged.
201  *
202  * Start by allocating an extra \ref GMX_SIMD_DOUBLE_WIDTH double elements of memory,
203  * and then call this function. The returned pointer will be greater or equal
204  * to the one you provided, and point to an address inside your provided memory
205  * that is aligned to the SIMD width.
206  */
207 static gmx_inline double *
208 gmx_simd_align_d(double *p)
209 {
210 #    ifdef GMX_SIMD_HAVE_DOUBLE
211     return (double *)(((size_t)((p)+GMX_SIMD_DOUBLE_WIDTH-1)) & (~((size_t)(GMX_SIMD_DOUBLE_WIDTH*sizeof(double)-1))));
212 #    else
213     return p;
214 #    endif
215 }
216
217 /*! \brief
218  * Align a (float) integer pointer for usage with SIMD instructions.
219  *
220  * You should typically \a not call this function directly (unless you explicitly
221  * want integers corresponding to single precision even when GMX_DOUBLE is
222  * set), but use the \ref gmx_simd_align_i macro to align integer memory
223  * corresponding to Gromacs default floating-point precision.
224  *
225  * \param  p Pointer to memory, allocate at least \ref GMX_SIMD_FINT32_WIDTH extra elements.
226  *
227  * \return Aligned pointer (>=p) suitable for loading/storing float-integer SIMD.
228  *         If \ref GMX_SIMD_HAVE_FINT32 is not set, p will be returned unchanged.
229  *
230  * This routine provides aligned memory for usage with \ref gmx_simd_fint32_t. You
231  * should have allocated an extra \ref GMX_SIMD_FINT32_WIDTH * sizeof(int) bytes. The
232  * reason why we need to separate float-integer vs. double-integer is that the
233  * width of registers after conversions from the floating-point types might not
234  * be identical, or even supported, in both cases.
235  */
236 static gmx_inline int *
237 gmx_simd_align_fi(int *p)
238 {
239 #    ifdef GMX_SIMD_HAVE_FINT32
240     return (int *)(((size_t)((p)+GMX_SIMD_FINT32_WIDTH-1)) & (~((size_t)(GMX_SIMD_FINT32_WIDTH*sizeof(int)-1))));
241 #    else
242     return p;
243 #    endif
244 }
245
246 /*! \brief
247  * Align a (double) integer pointer for usage with SIMD instructions.
248  *
249  * You should typically \a not call this function directly (unless you explicitly
250  * want integers corresponding to doublele precision even when GMX_DOUBLE is
251  * not set), but use the \ref gmx_simd_align_i macro to align integer memory
252  * corresponding to Gromacs default floating-point precision.
253  *
254  * \param  p Pointer to memory, allocate at least \ref GMX_SIMD_DINT32_WIDTH extra elements.
255  *
256  * \return Aligned pointer (>=p) suitable for loading/storing double-integer SIMD.
257  *         If \ref GMX_SIMD_HAVE_DINT32 is not set, p will be returned unchanged.
258  *
259  * This routine provides aligned memory for usage with \ref gmx_simd_dint32_t. You
260  * should have allocated an extra \ref GMX_SIMD_DINT32_WIDTH*sizeof(int) bytes. The
261  * reason why we need to separate float-integer vs. double-integer is that the
262  * width of registers after conversions from the floating-point types might not
263  * be identical, or even supported, in both cases.
264  */
265 static gmx_inline int *
266 gmx_simd_align_di(int *p)
267 {
268 #    ifdef GMX_SIMD_HAVE_DINT32
269     return (int *)(((size_t)((p)+GMX_SIMD_DINT32_WIDTH-1)) & (~((size_t)(GMX_SIMD_DINT32_WIDTH*sizeof(int)-1))));
270 #    else
271     return p;
272 #    endif
273 }
274
275 /*! \brief
276  * Align a float pointer for usage with SIMD4 instructions.
277  *
278  * You should typically \a not call this function directly (unless you explicitly
279  * want single precision even when GMX_DOUBLE is set), but use the
280  * \ref gmx_simd4_align_r macro to align memory in default Gromacs real precision.
281  *
282  * \param  p Pointer to memory, allocate at least \ref GMX_SIMD4_WIDTH extra elements.
283  *
284  * \return Aligned pointer (>=p) suitable for loading/storing float SIMD.
285  *         If \ref GMX_SIMD4_HAVE_FLOAT is not set, p will be returned unchanged.
286  *
287  * This routine provides aligned memory for usage with \ref gmx_simd4_float_t.
288  * should have allocated an extra \ref GMX_SIMD4_WIDTH * sizeof(float) bytes.
289  */
290 static gmx_inline float *
291 gmx_simd4_align_f(float *p)
292 {
293 #    ifdef GMX_SIMD4_HAVE_FLOAT
294     return (float *)(((size_t)((p)+GMX_SIMD4_WIDTH-1)) & (~((size_t)(GMX_SIMD4_WIDTH*sizeof(float)-1))));
295 #    else
296     return p;
297 #    endif
298 }
299
300 /*! \brief
301  * Align a double pointer for usage with SIMD4 instructions.
302  *
303  * You should typically \a not call this function directly (unless you explicitly
304  * want double precision even when GMX_DOUBLE is not set), but use the
305  * \ref gmx_simd4_align_r macro to align memory in default Gromacs real precision.
306  *
307  * \param  p Pointer to memory, allocate at least \ref GMX_SIMD4_WIDTH extra elements.
308  *
309  * \return Aligned pointer (>=p) suitable for loading/storing float SIMD.
310  *         If \ref GMX_SIMD4_HAVE_DOUBLE is not set, p will be returned unchanged.
311  *
312  * This routine provides aligned memory for usage with \ref gmx_simd4_double_t.
313  * should have allocated an extra \ref GMX_SIMD4_WIDTH * sizeof(double) bytes.
314  */
315 static gmx_inline double *
316 gmx_simd4_align_d(double *p)
317 {
318 #    ifdef GMX_SIMD4_HAVE_DOUBLE
319     return (double *)(((size_t)((p)+GMX_SIMD4_WIDTH-1)) & (~((size_t)(GMX_SIMD4_WIDTH*sizeof(double)-1))));
320 #    else
321     return p;
322 #    endif
323 }
324
325 /*! \} */
326
327
328 /* Define Gromacs "real" precision macros depending on Gromacs config. Note
329  * that conversions float-to-double and v.v. are not included here since they
330  * are not precision-dependent - find them in the implementation files.
331  */
332 #ifdef GMX_DOUBLE
333 /* Double floating-point. The documentation is in the float part below */
334 #    define gmx_simd_real_t                  gmx_simd_double_t
335 #    define gmx_simd_load_r                  gmx_simd_load_d
336 #    define gmx_simd_load1_r                 gmx_simd_load1_d
337 #    define gmx_simd_set1_r                  gmx_simd_set1_d
338 #    define gmx_simd_store_r                 gmx_simd_store_d
339 #    define gmx_simd_loadu_r                 gmx_simd_loadu_d
340 #    define gmx_simd_storeu_r                gmx_simd_storeu_d
341 #    define gmx_simd_setzero_r               gmx_simd_setzero_d
342 #    define gmx_simd_add_r                   gmx_simd_add_d
343 #    define gmx_simd_sub_r                   gmx_simd_sub_d
344 #    define gmx_simd_mul_r                   gmx_simd_mul_d
345 #    define gmx_simd_fmadd_r                 gmx_simd_fmadd_d
346 #    define gmx_simd_fmsub_r                 gmx_simd_fmsub_d
347 #    define gmx_simd_fnmadd_r                gmx_simd_fnmadd_d
348 #    define gmx_simd_fnmsub_r                gmx_simd_fnmsub_d
349 #    define gmx_simd_and_r                   gmx_simd_and_d
350 #    define gmx_simd_andnot_r                gmx_simd_andnot_d
351 #    define gmx_simd_or_r                    gmx_simd_or_d
352 #    define gmx_simd_xor_r                   gmx_simd_xor_d
353 #    define gmx_simd_rsqrt_r                 gmx_simd_rsqrt_d
354 #    define gmx_simd_rcp_r                   gmx_simd_rcp_d
355 #    define gmx_simd_fabs_r                  gmx_simd_fabs_d
356 #    define gmx_simd_fneg_r                  gmx_simd_fneg_d
357 #    define gmx_simd_max_r                   gmx_simd_max_d
358 #    define gmx_simd_min_r                   gmx_simd_min_d
359 #    define gmx_simd_round_r                 gmx_simd_round_d
360 #    define gmx_simd_trunc_r                 gmx_simd_trunc_d
361 #    define gmx_simd_fraction_r              gmx_simd_fraction_d
362 #    define gmx_simd_get_exponent_r          gmx_simd_get_exponent_d
363 #    define gmx_simd_get_mantissa_r          gmx_simd_get_mantissa_d
364 #    define gmx_simd_set_exponent_r          gmx_simd_set_exponent_d
365 /* Double integer and conversions */
366 #    define gmx_simd_int32_t                 gmx_simd_dint32_t
367 #    define gmx_simd_load_i                  gmx_simd_load_di
368 #    define gmx_simd_set1_i                  gmx_simd_set1_di
369 #    define gmx_simd_store_i                 gmx_simd_store_di
370 #    define gmx_simd_loadu_i                 gmx_simd_loadu_di
371 #    define gmx_simd_storeu_i                gmx_simd_storeu_di
372 #    define gmx_simd_setzero_i               gmx_simd_setzero_di
373 #    define gmx_simd_cvt_r2i                 gmx_simd_cvt_d2i
374 #    define gmx_simd_cvtt_r2i                gmx_simd_cvtt_d2i
375 #    define gmx_simd_cvt_i2r                 gmx_simd_cvt_i2d
376 #    define gmx_simd_extract_i               gmx_simd_extract_di
377 #    define gmx_simd_slli_i                  gmx_simd_slli_di
378 #    define gmx_simd_srli_i                  gmx_simd_srli_di
379 #    define gmx_simd_and_i                   gmx_simd_and_di
380 #    define gmx_simd_andnot_i                gmx_simd_andnot_di
381 #    define gmx_simd_or_i                    gmx_simd_or_di
382 #    define gmx_simd_xor_i                   gmx_simd_xor_di
383 #    define gmx_simd_add_i                   gmx_simd_add_di
384 #    define gmx_simd_sub_i                   gmx_simd_sub_di
385 #    define gmx_simd_mul_i                   gmx_simd_mul_di
386 /* Double booleans and selection */
387 #    define gmx_simd_bool_t                  gmx_simd_dbool_t
388 #    define gmx_simd_cmpeq_r                 gmx_simd_cmpeq_d
389 #    define gmx_simd_cmplt_r                 gmx_simd_cmplt_d
390 #    define gmx_simd_cmple_r                 gmx_simd_cmple_d
391 #    define gmx_simd_and_b                   gmx_simd_and_db
392 #    define gmx_simd_or_b                    gmx_simd_or_db
393 #    define gmx_simd_anytrue_b               gmx_simd_anytrue_db
394 #    define gmx_simd_blendzero_r             gmx_simd_blendzero_d
395 #    define gmx_simd_blendnotzero_r          gmx_simd_blendnotzero_d
396 #    define gmx_simd_blendv_r                gmx_simd_blendv_d
397 #    define gmx_simd_reduce_r                gmx_simd_reduce_d
398 #    define gmx_simd_ibool_t                 gmx_simd_dibool_t
399 #    define gmx_simd_cmpeq_i                 gmx_simd_cmpeq_di
400 #    define gmx_simd_cmplt_i                 gmx_simd_cmplt_di
401 #    define gmx_simd_and_ib                  gmx_simd_and_dib
402 #    define gmx_simd_or_ib                   gmx_simd_or_dib
403 #    define gmx_simd_anytrue_ib              gmx_simd_anytrue_dib
404 #    define gmx_simd_blendzero_i             gmx_simd_blendzero_di
405 #    define gmx_simd_blendnotzero_i          gmx_simd_blendnotzero_di
406 #    define gmx_simd_blendv_i                gmx_simd_blendv_di
407 /* Conversions between integer and double floating-point booleans */
408 #    define gmx_simd_cvt_b2ib                gmx_simd_cvt_db2dib
409 #    define gmx_simd_cvt_ib2b                gmx_simd_cvt_dib2db
410
411 /* SIMD4 double fp - we only support a subset of SIMD instructions for SIMD4 */
412 #    define gmx_simd4_real_t                 gmx_simd4_double_t
413 #    define gmx_simd4_load_r                 gmx_simd4_load_d
414 #    define gmx_simd4_load1_r                gmx_simd4_load1_d
415 #    define gmx_simd4_set1_r                 gmx_simd4_set1_d
416 #    define gmx_simd4_store_r                gmx_simd4_store_d
417 #    define gmx_simd4_loadu_r                gmx_simd4_loadu_d
418 #    define gmx_simd4_storeu_r               gmx_simd4_storeu_d
419 #    define gmx_simd4_setzero_r              gmx_simd4_setzero_d
420 #    define gmx_simd4_add_r                  gmx_simd4_add_d
421 #    define gmx_simd4_sub_r                  gmx_simd4_sub_d
422 #    define gmx_simd4_mul_r                  gmx_simd4_mul_d
423 #    define gmx_simd4_fmadd_r                gmx_simd4_fmadd_d
424 #    define gmx_simd4_fmsub_r                gmx_simd4_fmsub_d
425 #    define gmx_simd4_fnmadd_r               gmx_simd4_fnmadd_d
426 #    define gmx_simd4_fnmsub_r               gmx_simd4_fnmsub_d
427 #    define gmx_simd4_and_r                  gmx_simd4_and_d
428 #    define gmx_simd4_andnot_r               gmx_simd4_andnot_d
429 #    define gmx_simd4_or_r                   gmx_simd4_or_d
430 #    define gmx_simd4_xor_r                  gmx_simd4_xor_d
431 #    define gmx_simd4_rsqrt_r                gmx_simd4_rsqrt_d
432 #    define gmx_simd4_fabs_r                 gmx_simd4_fabs_d
433 #    define gmx_simd4_fneg_r                 gmx_simd4_fneg_d
434 #    define gmx_simd4_max_r                  gmx_simd4_max_d
435 #    define gmx_simd4_min_r                  gmx_simd4_min_d
436 #    define gmx_simd4_round_r                gmx_simd4_round_d
437 #    define gmx_simd4_trunc_r                gmx_simd4_trunc_d
438 #    define gmx_simd4_dotproduct3_r          gmx_simd4_dotproduct3_d
439 #    define gmx_simd4_bool_t                 gmx_simd4_dbool_t
440 #    define gmx_simd4_cmpeq_r                gmx_simd4_cmpeq_d
441 #    define gmx_simd4_cmplt_r                gmx_simd4_cmplt_d
442 #    define gmx_simd4_cmple_r                gmx_simd4_cmple_d
443 #    define gmx_simd4_and_b                  gmx_simd4_and_db
444 #    define gmx_simd4_or_b                   gmx_simd4_or_db
445 #    define gmx_simd4_anytrue_b              gmx_simd4_anytrue_db
446 #    define gmx_simd4_blendzero_r            gmx_simd4_blendzero_d
447 #    define gmx_simd4_blendnotzero_r         gmx_simd4_blendnotzero_d
448 #    define gmx_simd4_blendv_r               gmx_simd4_blendv_d
449 #    define gmx_simd4_reduce_r               gmx_simd4_reduce_d
450
451 /* Memory allocation */
452 #    define gmx_simd_align_r                 gmx_simd_align_d
453 #    define gmx_simd_align_i                 gmx_simd_align_di
454 #    define gmx_simd4_align_r                gmx_simd4_align_d
455
456 #    ifdef GMX_SIMD_HAVE_DOUBLE
457 #        define GMX_SIMD_HAVE_REAL
458 #        define GMX_SIMD_REAL_WIDTH          GMX_SIMD_DOUBLE_WIDTH
459 #    endif
460 #    ifdef GMX_SIMD_HAVE_DINT32
461 #        define GMX_SIMD_HAVE_INT32
462 #        define GMX_SIMD_INT32_WIDTH         GMX_SIMD_DINT32_WIDTH
463 #    endif
464 #    ifdef GMX_SIMD_HAVE_DINT32_EXTRACT
465 #        define GMX_SIMD_HAVE_INT32_EXTRACT
466 #    endif
467 #    ifdef GMX_SIMD_HAVE_DINT32_LOGICAL
468 #        define GMX_SIMD_HAVE_INT32_LOGICAL
469 #    endif
470 #    ifdef GMX_SIMD_HAVE_DINT32_ARITHMETICS
471 #        define GMX_SIMD_HAVE_INT32_ARITHMETICS
472 #    endif
473 #    ifdef GMX_SIMD4_HAVE_DOUBLE
474 #        define GMX_SIMD4_HAVE_REAL
475 #    endif
476
477 #else /* GMX_DOUBLE */
478
479 /*! \name SIMD data types
480  *
481  *  The actual storage of these types is implementation dependent. The
482  *  documentation is generated from the reference implementation, but for
483  *  normal usage this will likely not be what you are using.
484  * \{
485  */
486 /*! \brief Real precision floating-point SIMD datatype.
487  *
488  * This type is only available if \ref GMX_SIMD_HAVE_REAL is defined.
489  *
490  * If GMX_DOUBLE is defined, this will be set to \ref gmx_simd_double_t
491  * internally, otherwise \ref gmx_simd_float_t.
492  */
493 #    define gmx_simd_real_t                  gmx_simd_float_t
494
495 /*! \brief 32-bit integer SIMD type.
496  *
497  * This type is only available if \ref GMX_SIMD_HAVE_INT32 is defined.
498  *
499  * If GMX_DOUBLE is defined, this will be set to \ref gmx_simd_dint32_t
500  * internally, otherwise \ref gmx_simd_fint32_t. This might seem a strange
501  * implementation detail, but it is because some SIMD implementations use
502  * different types/widths of integers registers when converting from
503  * double vs. single precision floating point. As long as you just use
504  * this type you will not have to worry about precision.
505  */
506 #    define gmx_simd_int32_t                 gmx_simd_fint32_t
507
508 /*! \brief Boolean SIMD type for usage with \ref gmx_simd_real_t.
509  *
510  * This type is only available if \ref GMX_SIMD_HAVE_REAL is defined.
511  *
512  * If GMX_DOUBLE is defined, this will be set to \ref gmx_simd_dbool_t
513  * internally, otherwise \ref gmx_simd_fbool_t. This is necessary since some
514  * SIMD implementations use bitpatterns for marking truth, so single-
515  * vs. double precision booleans are not necessarily exchangable.
516  * As long as you just use this type you will not have to worry about precision.
517  *
518  * See \ref gmx_simd_ibool_t for an explanation of real vs. integer booleans.
519  */
520 #    define gmx_simd_bool_t                  gmx_simd_fbool_t
521
522 /*! \brief Boolean SIMD type for usage with \ref gmx_simd_int32_t.
523  *
524  * This type is only available if \ref GMX_SIMD_HAVE_INT32 is defined.
525  *
526  * If GMX_DOUBLE is defined, this will be set to \ref gmx_simd_dibool_t
527  * internally, otherwise \ref gmx_simd_fibool_t. This is necessary since some
528  * SIMD implementations use bitpatterns for marking truth, so single-
529  * vs. double precision booleans are not necessarily exchangable, and while
530  * a double-precision boolean might be represented with a 64-bit mask, the
531  * corresponding integer might only use a 32-bit mask.
532  *
533  * We provide conversion routines for these cases, so the only thing you need to
534  * keep in mind is to use \ref gmx_simd_bool_t when working with
535  * \ref gmx_simd_real_t while you pick \ref gmx_simd_ibool_t when working with
536  * \ref gmx_simd_int32_t.
537  *
538  * To convert between them, use \ref gmx_simd_cvt_b2ib and \ref gmx_simd_cvt_ib2b.
539  */
540 #    define gmx_simd_ibool_t                 gmx_simd_fibool_t
541
542
543 /*! \}
544  *  \name SIMD load/store operations on gmx_simd_real_t
545  *
546  *  \note Unaligned load/stores are only available when
547  *  \ref GMX_SIMD_HAVE_LOADU and \ref GMX_SIMD_HAVE_STOREU are set, respectively.
548  *  \{
549  */
550
551 /*! \brief Load \ref GMX_SIMD_REAL_WIDTH values from aligned memory to \ref gmx_simd_real_t
552  *
553  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_load_d,
554  * otherwise \ref gmx_simd_load_f.
555  *
556  * \copydetails gmx_simd_load_f
557  */
558 #    define gmx_simd_load_r                  gmx_simd_load_f
559
560 /*! \brief Set all elements in \ref gmx_simd_real_t from single value in memory.
561  *
562  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_load1_d,
563  * otherwise \ref gmx_simd_load1_f.
564  *
565  * \copydetails gmx_simd_load1_f
566  */
567 #    define gmx_simd_load1_r                 gmx_simd_load1_f
568
569 /*! \brief Set all elements in \ref gmx_simd_real_t from a scalar.
570  *
571  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_set1_d,
572  * otherwise \ref gmx_simd_set1_f.
573  *
574  * \copydetails gmx_simd_set1_f
575  */
576 #    define gmx_simd_set1_r                  gmx_simd_set1_f
577
578 /*! \brief Store \ref GMX_SIMD_REAL_WIDTH values from \ref gmx_simd_real_t to aligned memory.
579  *
580  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_store_d,
581  * otherwise \ref gmx_simd_store_f.
582  *
583  * \copydetails gmx_simd_store_f
584  */
585 #    define gmx_simd_store_r                 gmx_simd_store_f
586
587 /*! \brief Load \ref GMX_SIMD_REAL_WIDTH values from unaligned memory to \ref gmx_simd_real_t.
588  *
589  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_loadu_d,
590  * otherwise \ref gmx_simd_loadu_f.
591  *
592  * \copydetails gmx_simd_loadu_f
593  */
594 #    define gmx_simd_loadu_r                 gmx_simd_loadu_f
595
596 /*! \brief Store \ref GMX_SIMD_REAL_WIDTH values from \ref gmx_simd_real_t to unaligned memory.
597  *
598  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_storeu_d,
599  * otherwise \ref gmx_simd_storeu_f.
600  *
601  * \copydetails gmx_simd_storeu_f
602  */
603 #    define gmx_simd_storeu_r                gmx_simd_storeu_f
604
605 /*! \brief Set all elements in \ref gmx_simd_real_t to 0.0.
606  *
607  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_setzero_d,
608  * otherwise \ref gmx_simd_setzero_f.
609  *
610  * \copydetails gmx_simd_setzero_f
611  */
612 #    define gmx_simd_setzero_r               gmx_simd_setzero_f
613
614 /*! \}
615  *  \name SIMD load/store operations on gmx_simd_int32_t
616  *
617  *  \note Unaligned load/stores are only available when
618  *  \ref GMX_SIMD_HAVE_LOADU and \ref GMX_SIMD_HAVE_STOREU are set, respectively.
619  *  \{
620  */
621
622 /*! \brief Load \ref GMX_SIMD_INT32_WIDTH values from aligned memory to \ref gmx_simd_int32_t .
623  *
624  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_load_di ,
625  * otherwise \ref gmx_simd_load_fi .
626  *
627  * \copydetails gmx_simd_load_fi
628  */
629 #    define gmx_simd_load_i                  gmx_simd_load_fi
630
631 /*! \brief Set all elements in \ref gmx_simd_int32_t from a single integer.
632  *
633  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_set1_di ,
634  * otherwise \ref gmx_simd_set1_fi .
635  *
636  * \copydetails gmx_simd_set1_fi
637  */
638 #    define gmx_simd_set1_i                  gmx_simd_set1_fi
639
640 /*! \brief Store \ref GMX_SIMD_REAL_WIDTH values from \ref gmx_simd_int32_t to aligned memory.
641  *
642  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_store_di ,
643  * otherwise \ref gmx_simd_store_fi .
644  *
645  * \copydetails gmx_simd_store_fi
646  */
647 #    define gmx_simd_store_i                 gmx_simd_store_fi
648
649 /*! \brief Load \ref GMX_SIMD_REAL_WIDTH values from unaligned memory to \ref gmx_simd_int32_t.
650  *
651  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_loadu_di ,
652  * otherwise \ref gmx_simd_loadu_fi .
653  *
654  * \copydetails gmx_simd_loadu_fi
655  */
656 #    define gmx_simd_loadu_i                 gmx_simd_loadu_fi
657
658 /*! \brief Store \ref GMX_SIMD_REAL_WIDTH values from \ref gmx_simd_int32_t to unaligned memory.
659  *
660  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_storeu_di ,
661  * otherwise \ref gmx_simd_storeu_fi .
662  *
663  * \copydetails gmx_simd_storeu_fi
664  */
665 #    define gmx_simd_storeu_i                gmx_simd_storeu_fi
666
667 /*! \brief Extract single integer from \ref gmx_simd_int32_t element.
668  *
669  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_extract_di ,
670  * otherwise \ref gmx_simd_extract_fi .
671  *
672  * \copydetails gmx_simd_extract_fi
673  */
674 #    define gmx_simd_extract_i               gmx_simd_extract_fi
675
676 /*! \brief Set all elements in \ref gmx_simd_int32_t to 0.
677  *
678  * If GMX_DOUBLE is defined, it will be aliased to \ref gmx_simd_setzero_di ,
679  * otherwise \ref gmx_simd_setzero_fi .
680  *
681  * \copydetails gmx_simd_setzero_fi
682  */
683 #    define gmx_simd_setzero_i               gmx_simd_setzero_fi
684
685
686 /*! \}
687  *  \name SIMD floating-point logical operations on gmx_simd_real_t
688  *
689  *  These instructions are available if \ref GMX_SIMD_HAVE_LOGICAL is defined.
690  *  \{
691  */
692
693 /*! \brief Bitwise \a and on two \ref gmx_simd_real_t.
694  *
695  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_and_d,
696  * otherwise \ref gmx_simd_and_f.
697  *
698  * \copydetails gmx_simd_and_f
699  */
700 #    define gmx_simd_and_r                   gmx_simd_and_f
701
702 /*! \brief Bitwise \a and-not on two \ref gmx_simd_real_t; 1st arg is complemented.
703  *
704  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_andnot_d,
705  * otherwise \ref gmx_simd_andnot_f.
706  *
707  * \copydetails gmx_simd_andnot_f
708  */
709 #    define gmx_simd_andnot_r                gmx_simd_andnot_f
710
711 /*! \brief Bitwise \a or on two \ref gmx_simd_real_t.
712  *
713  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_or_d,
714  * otherwise \ref gmx_simd_or_f.
715  *
716  * \copydetails gmx_simd_or_f
717  */
718 #    define gmx_simd_or_r                    gmx_simd_or_f
719
720 /*! \brief Bitwise \a exclusive-or on two \ref gmx_simd_real_t.
721  *
722  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_xor_d,
723  * otherwise \ref gmx_simd_xor_f.
724  *
725  * \copydetails gmx_simd_xor_f
726  */
727 #    define gmx_simd_xor_r                   gmx_simd_xor_f
728
729 /*! \}
730  *  \name SIMD floating-point arithmetic operations on gmx_simd_real_t
731  *  \{
732  */
733
734 /*! \brief SIMD a+b for two \ref gmx_simd_real_t.
735  *
736  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_add_d,
737  * otherwise \ref gmx_simd_add_f.
738  *
739  * \copydetails gmx_simd_add_f
740  */
741 #    define gmx_simd_add_r                   gmx_simd_add_f
742
743 /*! \brief SIMD a-b for two \ref gmx_simd_real_t.
744  *
745  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_sub_d,
746  * otherwise \ref gmx_simd_sub_f.
747  *
748  * \copydetails gmx_simd_sub_f
749  */
750 #    define gmx_simd_sub_r                   gmx_simd_sub_f
751
752 /*! \brief SIMD a*b for two \ref gmx_simd_real_t.
753  *
754  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_mul_d,
755  * otherwise \ref gmx_simd_mul_f.
756  *
757  * \copydetails gmx_simd_mul_f
758  */
759 #    define gmx_simd_mul_r                   gmx_simd_mul_f
760
761 /*! \brief SIMD a*b+c for three \ref gmx_simd_real_t.
762  *
763  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_fmadd_d,
764  * otherwise \ref gmx_simd_fmadd_f.
765  *
766  * \copydetails gmx_simd_fmadd_f
767  */
768 #    define gmx_simd_fmadd_r                 gmx_simd_fmadd_f
769
770 /*! \brief SIMD a*b-c for three \ref gmx_simd_real_t.
771  *
772  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_fmsub_d,
773  * otherwise \ref gmx_simd_fmsub_f.
774  *
775  * \copydetails gmx_simd_fmsub_f
776  */
777 #    define gmx_simd_fmsub_r                 gmx_simd_fmsub_f
778
779 /*! \brief SIMD -a*b+c for three \ref gmx_simd_real_t.
780  *
781  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_fnmadd_d,
782  * otherwise \ref gmx_simd_fnmadd_f.
783  *
784  * \copydetails gmx_simd_fnmadd_f
785  */
786 #    define gmx_simd_fnmadd_r                gmx_simd_fnmadd_f
787
788 /*! \brief SIMD -a*b-c for three \ref gmx_simd_real_t.
789  *
790  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_fnmsub_d,
791  * otherwise \ref gmx_simd_fnmsub_f.
792  *
793  * \copydetails gmx_simd_fnmsub_f
794  */
795 #    define gmx_simd_fnmsub_r                gmx_simd_fnmsub_f
796
797 /*! \brief SIMD table lookup for 1/sqrt(x) approximation.
798  *
799  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_rsqrt_d,
800  * otherwise \ref gmx_simd_rsqrt_f.
801  *
802  * \copydetails gmx_simd_rsqrt_f
803  */
804 #    define gmx_simd_rsqrt_r                 gmx_simd_rsqrt_f
805
806 /*! \brief SIMD table lookup for 1/x approximation.
807  *
808  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_rcp_d,
809  * otherwise \ref gmx_simd_rcp_f.
810  *
811  * \copydetails gmx_simd_rcp_f
812  */
813 #    define gmx_simd_rcp_r                   gmx_simd_rcp_f
814
815 /*! \brief SIMD fabs(x) for \ref gmx_simd_real_t.
816  *
817  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_fabs_d,
818  * otherwise \ref gmx_simd_fabs_f.
819  *
820  * \copydetails gmx_simd_fabs_f
821  */
822 #    define gmx_simd_fabs_r                  gmx_simd_fabs_f
823
824 /*! \brief SIMD -x for \ref gmx_simd_real_t.
825  *
826  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_fneg_d,
827  * otherwise \ref gmx_simd_fneg_f.
828  *
829  * \copydetails gmx_simd_fneg_f
830  */
831 #    define gmx_simd_fneg_r                  gmx_simd_fneg_f
832
833 /*! \brief SIMD max(a,b) for each element in \ref gmx_simd_real_t.
834  *
835  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_max_d,
836  * otherwise \ref gmx_simd_max_f.
837  *
838  * \copydetails gmx_simd_max_f
839  */
840 #    define gmx_simd_max_r                   gmx_simd_max_f
841
842 /*! \brief SIMD min(a,b) for each element in \ref gmx_simd_real_t.
843  *
844  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_min_d,
845  * otherwise \ref gmx_simd_min_f.
846  *
847  * \copydetails gmx_simd_min_f
848  */
849 #    define gmx_simd_min_r                   gmx_simd_min_f
850
851 /*! \brief Round \ref gmx_simd_real_t to nearest int, return \ref gmx_simd_real_t.
852  *
853  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_round_d,
854  * otherwise \ref gmx_simd_round_f.
855  *
856  * \copydetails gmx_simd_round_f
857  */
858 #    define gmx_simd_round_r                 gmx_simd_round_f
859
860 /*! \brief Truncate \ref gmx_simd_real_t towards 0, return \ref gmx_simd_real_t.
861  *
862  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_trunc_d,
863  * otherwise \ref gmx_simd_trunc_f.
864  *
865  * \copydetails gmx_simd_trunc_f
866  */
867 #    define gmx_simd_trunc_r                 gmx_simd_trunc_f
868
869 /*! \brief SIMD Fraction, i.e. x-trunc(x) for \ref gmx_simd_real_t.
870  *
871  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_fraction_d,
872  * otherwise \ref gmx_simd_fraction_f.
873  *
874  * \copydetails gmx_simd_fraction_f
875  */
876 #    define gmx_simd_fraction_r              gmx_simd_fraction_f
877
878 /*! \brief Return the FP exponent of a SIMD \ref gmx_simd_real_t as a \ref gmx_simd_real_t.
879  *
880  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_get_exponent_d,
881  * otherwise \ref gmx_simd_get_exponent_f.
882  *
883  * \copydetails gmx_simd_exponent_f
884  */
885 #    define gmx_simd_get_exponent_r          gmx_simd_get_exponent_f
886
887 /*! \brief Return the FP mantissa of a SIMD \ref gmx_simd_real_t as a \ref gmx_simd_real_t.
888  *
889  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_get_mantissa_d,
890  * otherwise \ref gmx_simd_get_mantissa_f.
891  *
892  * \copydetails gmx_simd_mantissa_f
893  */
894 #    define gmx_simd_get_mantissa_r          gmx_simd_get_mantissa_f
895
896 /*! \brief Set the exponent of a SIMD \ref gmx_simd_real_t from a \ref gmx_simd_real_t.
897  *
898  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_set_exponent_d,
899  * otherwise \ref gmx_simd_set_exponent_f.
900  *
901  * \copydetails gmx_simd_set_exponent_f
902  */
903 #    define gmx_simd_set_exponent_r          gmx_simd_set_exponent_f
904
905 /*! \}
906  *  \name SIMD comparison, boolean, and select operations for gmx_simd_real_t
907  *  \{
908  */
909
910 /*! \brief SIMD a==b for \ref gmx_simd_real_t. Returns a \ref gmx_simd_bool_t.
911  *
912  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_cmpeq_d,
913  * otherwise \ref gmx_simd_cmpeq_f.
914  *
915  * \copydetails gmx_simd_cmpeq_f
916  */
917 #    define gmx_simd_cmpeq_r                 gmx_simd_cmpeq_f
918
919 /*! \brief SIMD a<b for \ref gmx_simd_real_t. Returns a \ref gmx_simd_bool_t.
920  *
921  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_cmplt_d,
922  * otherwise \ref gmx_simd_cmplt_f.
923  *
924  * \copydetails gmx_simd_cmplt_f
925  */
926 #    define gmx_simd_cmplt_r                 gmx_simd_cmplt_f
927
928 /*! \brief SIMD a<=b for \ref gmx_simd_real_t. Returns a \ref gmx_simd_bool_t.
929  *
930  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_cmple_d,
931  * otherwise \ref gmx_simd_cmple_f.
932  *
933  * \copydetails gmx_simd_cmple_f
934  */
935 #    define gmx_simd_cmple_r                 gmx_simd_cmple_f
936
937 /*! \brief For each element, the result boolean is true if both arguments are true
938  *
939  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_and_db,
940  * otherwise \ref gmx_simd_and_fb.
941  *
942  * \copydetails gmx_simd_and_fb
943  */
944 #    define gmx_simd_and_b                   gmx_simd_and_fb
945
946 /*! \brief For each element, the result boolean is true if either argument is true
947  *
948  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_or_db,
949  * otherwise \ref gmx_simd_or_fb.
950  *
951  * \copydetails gmx_simd_or_fn
952  */
953 #    define gmx_simd_or_b                    gmx_simd_or_fb
954
955 /*! \brief Return nonzero if any element in gmx_simd_bool_t is true, otherwise 0.
956  *
957  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_anytrue_db,
958  * otherwise \ref gmx_simd_anytrue_fb.
959  *
960  * \copydetails gmx_simd_anytrue_fb
961  */
962 #    define gmx_simd_anytrue_b               gmx_simd_anytrue_fb
963
964 /*! \brief Selects elements from \ref gmx_simd_real_t where boolean is true, otherwise 0.
965  *
966  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_blendzero_d,
967  * otherwise \ref gmx_simd_blendzero_f.
968  *
969  * \copydetails gmx_simd_blendzero_f
970  *
971  * \sa gmx_simd_blendzero_i
972  */
973 #    define gmx_simd_blendzero_r             gmx_simd_blendzero_f
974
975 /*! \brief Selects elements from \ref gmx_simd_real_t where boolean is false, otherwise 0.
976  *
977  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_blendnotzero_d,
978  * otherwise \ref gmx_simd_blendnotzero_f.
979  *
980  * \copydetails gmx_simd_blendnotzero_f
981  */
982 #    define gmx_simd_blendnotzero_r          gmx_simd_blendnotzero_f
983
984 /*! \brief Selects from 2nd real SIMD arg where boolean is true, otherwise 1st arg.
985  *
986  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_blendv_d,
987  * otherwise \ref gmx_simd_blendv_f.
988  *
989  * \copydetails gmx_simd_blendv_f
990  */
991 #    define gmx_simd_blendv_r                gmx_simd_blendv_f
992
993 /*! \brief Return sum of all elements in SIMD floating-point variable.
994  *
995  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_reduce_d,
996  * otherwise \ref gmx_simd_reduce_f.
997  *
998  * \copydetails gmx_simd_reduce_f
999  */
1000 #    define gmx_simd_reduce_r                gmx_simd_reduce_f
1001
1002 /*! \}
1003  *  \name SIMD integer logical operations on gmx_simd_int32_t
1004  *
1005  *  These instructions are available if \ref GMX_SIMD_HAVE_INT32_LOGICAL is defined.
1006  *  \{
1007  */
1008
1009 /*! \brief Shift each element in \ref gmx_simd_int32_t left by immediate
1010  *
1011  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_slli_di,
1012  * otherwise \ref gmx_simd_slli_fi.
1013  *
1014  * \copydetails gmx_simd_slli_fi
1015  */
1016 #    define gmx_simd_slli_i                  gmx_simd_slli_fi
1017
1018 /*! \brief Shift each element in \ref gmx_simd_int32_t right by immediate
1019  *
1020  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_srli_di,
1021  * otherwise \ref gmx_simd_srli_fi.
1022  *
1023  * \copydetails gmx_simd_srli_fi
1024  */
1025 #    define gmx_simd_srli_i                  gmx_simd_srli_fi
1026
1027 /*! \brief Bitwise \a and on two \ref gmx_simd_int32_t.
1028  *
1029  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_and_di,
1030  * otherwise \ref gmx_simd_and_fi.
1031  *
1032  * \copydetails gmx_simd_and_fi
1033  */
1034 #    define gmx_simd_and_i                   gmx_simd_and_fi
1035
1036 /*! \brief Bitwise \a and-not on two \ref gmx_simd_int32_t; 1st arg is complemented.
1037  *
1038  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_andnot_di,
1039  * otherwise \ref gmx_simd_andnot_fi.
1040  *
1041  * \copydetails gmx_simd_andnot_fi
1042  */
1043 #    define gmx_simd_andnot_i                gmx_simd_andnot_fi
1044
1045 /*! \brief Bitwise \a or on two \ref gmx_simd_int32_t.
1046  *
1047  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_or_di,
1048  * otherwise \ref gmx_simd_or_fi.
1049  *
1050  * \copydetails gmx_simd_or_fi
1051  */
1052 #    define gmx_simd_or_i                    gmx_simd_or_fi
1053
1054 /*! \brief Bitwise \a xor on two \ref gmx_simd_int32_t.
1055  *
1056  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_xor_di,
1057  * otherwise \ref gmx_simd_xor_fi.
1058  *
1059  * \copydetails gmx_simd_xor_fi
1060  */
1061 #    define gmx_simd_xor_i                   gmx_simd_xor_fi
1062
1063 /*! \}
1064  *  \name SIMD integer arithmetic operations on gmx_simd_int32_t
1065  *
1066  *  These instructions are available if \ref GMX_SIMD_HAVE_INT32_ARITHMETICS is defined.
1067  *  \{
1068  */
1069
1070 /*! \brief SIMD a+b for two \ref gmx_simd_int32_t.
1071  *
1072  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_add_di,
1073  * otherwise \ref gmx_simd_add_fi.
1074  *
1075  * \copydetails gmx_simd_add_fi
1076  */
1077 #    define gmx_simd_add_i                   gmx_simd_add_fi
1078
1079 /*! \brief SIMD a-b for two \ref gmx_simd_int32_t.
1080  *
1081  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_sub_di,
1082  * otherwise \ref gmx_simd_sub_fi.
1083  *
1084  * \copydetails gmx_simd_sub_fi
1085  */
1086 #    define gmx_simd_sub_i                   gmx_simd_sub_fi
1087
1088 /*! \brief SIMD a*b for two \ref gmx_simd_int32_t.
1089  *
1090  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_mul_di,
1091  * otherwise \ref gmx_simd_mul_fi.
1092  *
1093  * \copydetails gmx_simd_mul_fi
1094  */
1095 #    define gmx_simd_mul_i                   gmx_simd_mul_fi
1096
1097 /*! \}
1098  *  \name SIMD integer comparison, booleans, and selection on gmx_simd_int32_t
1099  *
1100  *  These instructions are available if \ref GMX_SIMD_HAVE_INT32_ARITHMETICS is defined.
1101  *  \{
1102  */
1103
1104 /*! \brief Returns boolean describing whether a==b, for \ref gmx_simd_int32_t
1105  *
1106  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_cmpeq_di,
1107  * otherwise \ref gmx_simd_cmpeq_fi.
1108  *
1109  * \copydetails gmx_simd_cmpeq_fi
1110  */
1111 #    define gmx_simd_cmpeq_i                 gmx_simd_cmpeq_fi
1112
1113 /*! \brief Returns boolean describing whether a<b, for \ref gmx_simd_int32_t
1114  *
1115  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_cmplt_di,
1116  * otherwise \ref gmx_simd_cmplt_fi.
1117  *
1118  * \copydetails gmx_simd_cmplt_fi
1119  */
1120 #    define gmx_simd_cmplt_i                 gmx_simd_cmplt_fi
1121
1122 /*! \brief For each element, the result boolean is true if both arguments are true
1123  *
1124  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_and_dib,
1125  * otherwise \ref gmx_simd_and_fib.
1126  *
1127  * \copydetails gmx_simd_and_fib
1128  */
1129 #    define gmx_simd_and_ib                  gmx_simd_and_fib
1130
1131 /*! \brief For each element, the result boolean is true if either argument is true.
1132  *
1133  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_or_dib,
1134  * otherwise \ref gmx_simd_or_fib.
1135  *
1136  * \copydetails gmx_simd_or_fib
1137  */
1138 #    define gmx_simd_or_ib                   gmx_simd_or_fib
1139
1140 /*! \brief Return nonzero if any element in gmx_simd_ibool_t is true, otherwise 0.
1141  *
1142  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_anytrue_dib,
1143  * otherwise \ref gmx_simd_anytrue_fib.
1144  *
1145  * \copydetails gmx_simd_anytrue_fib
1146  */
1147 #    define gmx_simd_anytrue_ib              gmx_simd_anytrue_fib
1148
1149 /*! \brief Selects elements from \ref gmx_simd_int32_t where boolean is true, otherwise 0.
1150  *
1151  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_blendzero_di,
1152  * otherwise \ref gmx_simd_blendzero_fi.
1153  *
1154  * \copydetails gmx_simd_blendzero_fi
1155  */
1156 #    define gmx_simd_blendzero_i             gmx_simd_blendzero_fi
1157
1158 /*! \brief Selects elements from \ref gmx_simd_int32_t where boolean is false, otherwise 0.
1159  *
1160  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_blendnotzero_di,
1161  * otherwise \ref gmx_simd_blendnotzero_fi.
1162  *
1163  * \copydetails gmx_simd_blendnotzero_fi
1164  */
1165 #    define gmx_simd_blendnotzero_i          gmx_simd_blendnotzero_fi
1166
1167 /*! \brief Selects from 2nd int SIMD arg where boolean is true, otherwise 1st arg.
1168  *
1169  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_blendv_di,
1170  * otherwise \ref gmx_simd_blendv_fi.
1171  *
1172  * \copydetails gmx_simd_blendv_fi
1173  */
1174 #    define gmx_simd_blendv_i                gmx_simd_blendv_fi
1175
1176 /*! \}
1177  *  \name SIMD conversion operations
1178  *
1179  *  These instructions are available when both types involved in the conversion
1180  *  are defined, e.g. \ref GMX_SIMD_HAVE_REAL and \ref GMX_SIMD_HAVE_INT32
1181  *  for real-to-integer conversion.
1182  *  \{
1183  */
1184
1185 /*! \brief Convert gmx_simd_real_t to gmx_simd_int32_t, round to nearest integer.
1186  *
1187  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_cvt_d2i,
1188  * otherwise \ref gmx_simd_cvt_f2i.
1189  *
1190  * \copydetails gmx_simd_cvt_f2i
1191  */
1192 #    define gmx_simd_cvt_r2i                 gmx_simd_cvt_f2i
1193
1194 /*! \brief Convert gmx_simd_real_t to gmx_simd_int32_t, truncate towards zero
1195  *
1196  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_cvtt_d2i,
1197  * otherwise \ref gmx_simd_cvtt_f2i.
1198  *
1199  * \copydetails gmx_simd_cvtt_f2i
1200  */
1201 #    define gmx_simd_cvtt_r2i                gmx_simd_cvtt_f2i
1202
1203 /*! \brief Convert gmx_simd_int32_t to gmx_simd_real_t
1204  *
1205  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_cvt_i2d,
1206  * otherwise \ref gmx_simd_cvt_i2f.
1207  *
1208  * \copydetails gmx_simd_cvt_i2f
1209  */
1210 #    define gmx_simd_cvt_i2r                 gmx_simd_cvt_i2f
1211
1212 /*! \brief Convert from gmx_simd_bool_t to gmx_simd_ibool_t
1213  *
1214  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_cvt_db2dib,
1215  * otherwise \ref gmx_simd_cvt_fb2fib.
1216  *
1217  * \copydetails gmx_simd_cvt_fb2fib
1218  */
1219 #    define gmx_simd_cvt_b2ib                gmx_simd_cvt_fb2fib
1220
1221 /*! \brief Convert from gmx_simd_ibool_t to gmx_simd_bool_t
1222  *
1223  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_cvt_dib2db,
1224  * otherwise \ref gmx_simd_cvt_fib2fb.
1225  *
1226  * \copydetails gmx_simd_cvt_fib2fb
1227  */
1228 #    define gmx_simd_cvt_ib2b                gmx_simd_cvt_fib2fb
1229
1230
1231 /*! \}
1232  *  \name SIMD memory alignment operations
1233  *  \{
1234  */
1235
1236 /*! \brief Align real memory for SIMD usage.
1237  *
1238  * This routine will only align memory if \ref GMX_SIMD_HAVE_REAL is defined.
1239  * Otherwise the original pointer will be returned.
1240  *
1241  * Start by allocating an extra \ref GMX_SIMD_REAL_WIDTH float elements of memory,
1242  * and then call this function. The returned pointer will be greater or equal
1243  * to the one you provided, and point to an address inside your provided memory
1244  * that is aligned to the SIMD width.
1245  *
1246  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_align_d,
1247  * otherwise \ref gmx_simd_align_f. For detailed documentation, see the
1248  * precision-specific implementation routines.
1249  */
1250 #    define gmx_simd_align_r                 gmx_simd_align_f
1251
1252 /*! \brief Align integer memory for SIMD usage.
1253  *
1254  * This routine will only align memory if \ref GMX_SIMD_HAVE_INT32 is defined.
1255  * Otherwise the original pointer will be returned.
1256  *
1257  * Start by allocating an extra \ref GMX_SIMD_INT32_WIDTH elements of memory,
1258  * and then call this function. The returned pointer will be greater or equal
1259  * to the one you provided, and point to an address inside your provided memory
1260  * that is aligned to the SIMD width.
1261  *
1262  * If GMX_DOUBLE is defined, this will be aliased to \ref gmx_simd_align_di,
1263  * otherwise \ref gmx_simd_align_fi. For detailed documentation, see the
1264  * precision-specific implementation routines.
1265  */
1266 #    define gmx_simd_align_i                 gmx_simd_align_fi
1267
1268 /*! \} */
1269
1270 /*! \name SIMD4 - constant width-four SIMD datatypes
1271  *
1272  * These operations are only meant to be used for a few coordinate
1273  * manipulation and grid interpolation routines, so we only support a subset
1274  * of operations for SIMD4. To avoid repeating all the documentation from
1275  * the generic width SIMD routines, we only provide brief documentation for
1276  * these operations. Follow the link to the implementation documentation or the
1277  * reference to the corresponding generic SIMD routine. The format will be
1278  * exactly the same, but they have SIMD replaced with SIMD4.
1279  *  \{
1280  */
1281
1282 /*! \brief SIMD real datatype guaranteed to be 4 elements wide, if available.
1283  *
1284  * All the SIMD4 datatypes and operations behave like their counterparts for
1285  * the generic SIMD implementation, but they might be implemented with different
1286  * registers, or not supported at all. It is important that you check the
1287  * define \ref GMX_SIMD4_HAVE_REAL before using it.
1288  *
1289  * Just as the normal SIMD operations, all SIMD4 types and routines will
1290  * be aliased to either single or double precision ones based on whether
1291  * GMX_DOUBLE is defined.
1292  *
1293  * \note There is no support for integer or math operations in SIMD4.
1294  */
1295 #    define gmx_simd4_real_t                 gmx_simd4_float_t
1296
1297 /*! \brief Boolean for \ref gmx_simd4_real_t comparision/selection */
1298 #    define gmx_simd4_bool_t                 gmx_simd4_fbool_t
1299
1300 /*! \brief Load aligned data to gmx_simd4_real_t.
1301  *
1302  * \copydetails gmx_simd4_load_f
1303  */
1304 #    define gmx_simd4_load_r                 gmx_simd4_load_f
1305
1306 /*! \brief Load single element to gmx_simd4_real_t
1307  *
1308  * \copydetails gmx_simd4_load1_f
1309  */
1310 #    define gmx_simd4_load1_r                gmx_simd4_load1_f
1311
1312 /*! \brief Set gmx_simd4_real_t from scalar value
1313  *
1314  * \copydetails gmx_simd4_set1_f
1315  */
1316 #    define gmx_simd4_set1_r                 gmx_simd4_set1_f
1317
1318 /*! \brief store aligned data from gmx_simd4_real_t
1319  *
1320  * \copydetails gmx_simd4_store_f
1321  */
1322 #    define gmx_simd4_store_r                gmx_simd4_store_f
1323
1324 /*! \brief Load unaligned data to gmx_simd4_real_t
1325  *
1326  * \copydetails gmx_simd4_loadu_f
1327  */
1328 #    define gmx_simd4_loadu_r                gmx_simd4_loadu_f
1329
1330 /*! \brief Store unaligned data from gmx_simd4_real_t
1331  *
1332  * \copydetails gmx_simd4_storeu_f
1333  */
1334 #    define gmx_simd4_storeu_r               gmx_simd4_storeu_f
1335
1336 /*! \brief Set all elements in gmx_simd4_real_t to 0.0
1337  *
1338  * \copydetails gmx_simd4_setzero_f
1339  */
1340 #    define gmx_simd4_setzero_r              gmx_simd4_setzero_f
1341
1342 /*! \brief Bitwise and for two gmx_simd4_real_t
1343  *
1344  * \copydetails gmx_simd4_and_f
1345  */
1346 #    define gmx_simd4_and_r                  gmx_simd4_and_f
1347
1348 /*! \brief Bitwise and-not for two gmx_simd4_real_t. 1st arg is complemented.
1349  *
1350  * \copydetails gmx_simd4_andnot_f
1351  */
1352 #    define gmx_simd4_andnot_r               gmx_simd4_andnot_f
1353
1354 /*! \brief Bitwise or for two gmx_simd4_real_t
1355  *
1356  * \copydetails gmx_simd4_or_f
1357  */
1358 #    define gmx_simd4_or_r                   gmx_simd4_or_f
1359
1360 /*! \brief Bitwise xor for two gmx_simd4_real_t
1361  *
1362  * \copydetails gmx_simd4_xor_f
1363  */
1364 #    define gmx_simd4_xor_r                  gmx_simd4_xor_f
1365
1366 /*! \brief a+b for \ref gmx_simd4_real_t
1367  *
1368  * \copydetails gmx_simd4_add_f
1369  */
1370 #    define gmx_simd4_add_r                  gmx_simd4_add_f
1371
1372 /*! \brief a-b for \ref gmx_simd4_real_t
1373  *
1374  * \copydetails gmx_simd4_sub_f
1375  */
1376 #    define gmx_simd4_sub_r                  gmx_simd4_sub_f
1377
1378 /*! \brief a*b for \ref gmx_simd4_real_t
1379  *
1380  * \copydetails gmx_simd4_mul_f
1381  */
1382 #    define gmx_simd4_mul_r                  gmx_simd4_mul_f
1383
1384 /*! \brief a*b+c for \ref gmx_simd4_real_t
1385  *
1386  * \copydetails gmx_simd4_fmadd_f
1387  */
1388 #    define gmx_simd4_fmadd_r                gmx_simd4_fmadd_f
1389
1390 /*! \brief a*b-c for \ref gmx_simd4_real_t
1391  *
1392  * \copydetails gmx_simd4_fmsub_f
1393  */
1394 #    define gmx_simd4_fmsub_r                gmx_simd4_fmsub_f
1395
1396 /*! \brief -a*b+c for \ref gmx_simd4_real_t
1397  *
1398  * \copydetails gmx_simd4_fnmadd_f
1399  */
1400 #    define gmx_simd4_fnmadd_r               gmx_simd4_fnmadd_f
1401
1402 /*! \brief -a*b-c for \ref gmx_simd4_real_t
1403  *
1404  * \copydetails gmx_simd4_fnmsub_f
1405  */
1406 #    define gmx_simd4_fnmsub_r               gmx_simd4_fnmsub_f
1407
1408 /*! \brief 1/sqrt(x) approximate lookup for \ref gmx_simd4_real_t
1409  *
1410  * \copydetails gmx_simd4_rsqrt_f
1411  */
1412 #    define gmx_simd4_rsqrt_r                gmx_simd4_rsqrt_f
1413
1414 /*! \brief fabs(x) for \ref gmx_simd4_real_t
1415  *
1416  * \copydetails gmx_simd4_fabs_f
1417  */
1418 #    define gmx_simd4_fabs_r                 gmx_simd4_fabs_f
1419
1420 /*! \brief Change sign (-x) for \ref gmx_simd4_real_t
1421  *
1422  * \copydetails gmx_simd4_fneg_f
1423  */
1424 #    define gmx_simd4_fneg_r                 gmx_simd4_fneg_f
1425
1426 /*! \brief Select maximum of each pair of elements from args for \ref gmx_simd4_real_t
1427  *
1428  * \copydetails gmx_simd4_max_f
1429  */
1430 #    define gmx_simd4_max_r                  gmx_simd4_max_f
1431
1432 /*! \brief Select minimum of each pair of elements from args for \ref gmx_simd4_real_t
1433  *
1434  * \copydetails gmx_simd4_min_f
1435  */
1436 #    define gmx_simd4_min_r                  gmx_simd4_min_f
1437
1438 /*! \brief Round \ref gmx_simd4_real_t to nearest integer, return \ref gmx_simd4_real_t
1439  *
1440  * \copydetails gmx_simd4_round_f
1441  */
1442 #    define gmx_simd4_round_r                gmx_simd4_round_f
1443
1444 /*! \brief Truncate \ref gmx_simd4_real_t towards zero, return \ref gmx_simd4_real_t
1445  *
1446  * \copydetails gmx_simd4_trunc_f
1447  */
1448 #    define gmx_simd4_trunc_r                gmx_simd4_trunc_f
1449
1450 /*! \brief Scalar product of first three elements of two \ref gmx_simd4_real_t *
1451  *
1452  * \copydetails gmx_simd4_dotproduct3_f
1453  */
1454 #    define gmx_simd4_dotproduct3_r          gmx_simd4_dotproduct3_f
1455
1456 /*! \brief Return booleans whether a==b for each element two \ref gmx_simd4_real_t
1457  *
1458  * \copydetails gmx_simd4_cmpeq_f
1459  */
1460 #    define gmx_simd4_cmpeq_r                gmx_simd4_cmpeq_f
1461 /*! \brief Return booleans whether a<b for each element two \ref gmx_simd4_real_t
1462  *
1463  * \copydetails gmx_simd4_cmplt_f
1464  */
1465 #    define gmx_simd4_cmplt_r                gmx_simd4_cmplt_f
1466 /*! \brief Return booleans whether a<=b for each element two \ref gmx_simd4_real_t
1467  *
1468  * \copydetails gmx_simd4_cmple_f
1469  */
1470 #    define gmx_simd4_cmple_r                gmx_simd4_cmple_f
1471
1472 /*! \brief Logical and for two \ref gmx_simd4_bool_t
1473  *
1474  * \copydetails gmx_simd4_and_fb
1475  */
1476 #    define gmx_simd4_and_b                  gmx_simd4_and_fb
1477 /*! \brief Logical or for two \ref gmx_simd4_bool_t
1478  *
1479  * \copydetails gmx_simd4_or_fb
1480  */
1481 #    define gmx_simd4_or_b                   gmx_simd4_or_fb
1482
1483 /*! \brief Return nonzero if any element in \ref gmx_simd4_bool_t is true, otherwise 0
1484  *
1485  * \copydetails gmx_simd4_anytrue_fb
1486  */
1487 #    define gmx_simd4_anytrue_b              gmx_simd4_anytrue_fb
1488
1489 /*! \brief Selects from 2nd real SIMD4 arg where boolean is true, otherwise 1st arg
1490  *
1491  * \copydetails gmx_simd4_blendzero_f
1492  */
1493 #    define gmx_simd4_blendzero_r            gmx_simd4_blendzero_f
1494
1495 /*! \brief Selects from 2nd real SIMD4 arg where boolean is false, otherwise 1st arg
1496  *
1497  * \copydetails gmx_simd4_blendnotzero_f
1498  */
1499 #    define gmx_simd4_blendnotzero_r            gmx_simd4_blendnotzero_f
1500
1501 /*! \brief Selects from 2nd real SIMD4 arg where boolean is true, otherwise 1st arg
1502  *
1503  * \copydetails gmx_simd4_blendv_f
1504  */
1505 #    define gmx_simd4_blendv_r               gmx_simd4_blendv_f
1506
1507 /*! \brief Return sum of all elements in SIMD4 floating-point variable.
1508  *
1509  * \copydetails gmx_simd4_reduce_f
1510  */
1511 #    define gmx_simd4_reduce_r               gmx_simd4_reduce_f
1512
1513 /*! \brief Align real memory for SIMD4 usage.
1514  *
1515  * \copydetails gmx_simd4_align_f
1516  */
1517 #    define gmx_simd4_align_r                gmx_simd4_align_f
1518
1519 /*! \} */
1520
1521 /*! \name SIMD predefined macros to describe high-level capabilities
1522  *  \{
1523  */
1524
1525 #    if (defined GMX_SIMD_HAVE_FLOAT) || (defined DOXYGEN)
1526 /*! \brief Defined if gmx_simd_real_t is available.
1527  *
1528  *  if GMX_DOUBLE is defined, this will be aliased to
1529  *  \ref GMX_SIMD_HAVE_DOUBLE, otherwise GMX_SIMD_HAVE_FLOAT.
1530  */
1531 #        define GMX_SIMD_HAVE_REAL
1532 /*! \brief Width of gmx_simd_real_t.
1533  *
1534  *  if GMX_DOUBLE is defined, this will be aliased to
1535  *  \ref GMX_SIMD_DOUBLE_WIDTH, otherwise GMX_SIMD_FLOAT_WIDTH.
1536  */
1537 #        define GMX_SIMD_REAL_WIDTH          GMX_SIMD_FLOAT_WIDTH
1538 #    endif
1539 #    if (defined GMX_SIMD_HAVE_FINT32) || (defined DOXYGEN)
1540 /*! \brief Defined if gmx_simd_int32_t is available.
1541  *
1542  *  if GMX_DOUBLE is defined, this will be aliased to
1543  *  \ref GMX_SIMD_HAVE_DINT32, otherwise GMX_SIMD_HAVE_FINT32.
1544  */
1545 #        define GMX_SIMD_HAVE_INT32
1546 /*! \brief Width of gmx_simd_int32_t.
1547  *
1548  *  if GMX_DOUBLE is defined, this will be aliased to
1549  *  \ref GMX_SIMD_DINT32_WIDTH, otherwise GMX_SIMD_FINT32_WIDTH.
1550  */
1551 #        define GMX_SIMD_INT32_WIDTH         GMX_SIMD_FINT32_WIDTH
1552 #    endif
1553 #    if (defined GMX_SIMD_HAVE_FINT32_EXTRACT) || (defined DOXYGEN)
1554 /*! \brief Defined if gmx_simd_extract_i() is available.
1555  *
1556  *  if GMX_DOUBLE is defined, this will be aliased to
1557  *  \ref GMX_SIMD_HAVE_DINT32_EXTRACT, otherwise GMX_SIMD_HAVE_FINT32_EXTRACT.
1558  */
1559 #        define GMX_SIMD_HAVE_INT32_EXTRACT
1560 #    endif
1561 #    if (defined GMX_SIMD_HAVE_FINT32_LOGICAL) || (defined DOXYGEN)
1562 /*! \brief Defined if logical ops are supported on gmx_simd_int32_t.
1563  *
1564  *  if GMX_DOUBLE is defined, this will be aliased to
1565  *  \ref GMX_SIMD_HAVE_DINT32_LOGICAL, otherwise GMX_SIMD_HAVE_FINT32_LOGICAL.
1566  */
1567 #        define GMX_SIMD_HAVE_INT32_LOGICAL
1568 #    endif
1569 #    if (defined GMX_SIMD_HAVE_FINT32_ARITHMETICS) || (defined DOXYGEN)
1570 /*! \brief Defined if arithmetic ops are supported on gmx_simd_int32_t.
1571  *
1572  *  if GMX_DOUBLE is defined, this will be aliased to
1573  *  \ref GMX_SIMD_HAVE_DINT32_ARITHMETICS, otherwise GMX_SIMD_HAVE_FINT32_ARITHMETICS.
1574  */
1575 #        define GMX_SIMD_HAVE_INT32_ARITHMETICS
1576 #    endif
1577 #    if (defined GMX_SIMD4_HAVE_FLOAT) || (defined DOXYGEN)
1578 /*! \brief Defined if gmx_simd4_real_t is available.
1579  *
1580  *  if GMX_DOUBLE is defined, this will be aliased to
1581  *  \ref GMX_SIMD4_HAVE_DOUBLE, otherwise GMX_SIMD4_HAVE_FLOAT.
1582  */
1583 #        define GMX_SIMD4_HAVE_REAL
1584 #    endif
1585
1586 /*! \} */
1587
1588 #endif /* GMX_DOUBLE */
1589
1590 /*! \} */
1591 /*! \endcond */
1592
1593 #endif /* GMX_SIMD_SIMD_H */