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