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