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