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