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