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