First part of commit for redesigned SIMD module - namechanges.
[alexxy/gromacs.git] / src / gromacs / simd / macros_ref.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 #ifndef GMX_SIMD_MACROS_REF_H
37 #define GMX_SIMD_MACROS_REF_H
38
39 /* This file contains a reference plain-C implementation of arbitrary width.
40  * This code is only useful for testing and documentation.
41  * The SIMD width is set by defining GMX_SIMD_REF_WIDTH before including.
42  */
43
44
45 #ifndef GMX_SIMD_REF_WIDTH
46 #error "GMX_SIMD_REF_WIDTH should be defined before including gromacs/simd/macros_ref.h"
47 #endif
48
49 #include <math.h>
50
51 /* float/double SIMD register type */
52 typedef struct {
53     real r[GMX_SIMD_REF_WIDTH];
54 } gmx_simd_ref_pr;
55
56 /* boolean SIMD register type */
57 typedef struct {
58     char r[GMX_SIMD_REF_WIDTH];
59 } gmx_simd_ref_pb;
60
61 /* integer SIMD register type, only for table indexing and exclusion masks */
62 typedef struct {
63     int r[GMX_SIMD_REF_WIDTH];
64 } gmx_simd_ref_epi32;
65 #define GMX_SIMD_REF_EPI32_WIDTH  GMX_SIMD_REF_WIDTH
66
67 /* Load GMX_SIMD_REF_WIDTH reals for memory starting at r */
68 static gmx_inline gmx_simd_ref_pr
69 gmx_simd_ref_load_pr(const real *r)
70 {
71     gmx_simd_ref_pr a;
72     int             i;
73
74     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
75     {
76         a.r[i] = r[i];
77     }
78
79     return a;
80 }
81
82 /* Set all SIMD register elements to *r */
83 static gmx_inline gmx_simd_ref_pr
84 gmx_simd_ref_load1_pr(const real *r)
85 {
86     gmx_simd_ref_pr a;
87     int             i;
88
89     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
90     {
91         a.r[i] = *r;
92     }
93
94     return a;
95 }
96
97 /* Set all SIMD register elements to r */
98 static gmx_inline gmx_simd_ref_pr
99 gmx_simd_ref_set1_pr(real r)
100 {
101     gmx_simd_ref_pr a;
102     int             i;
103
104     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
105     {
106         a.r[i] = r;
107     }
108
109     return a;
110 }
111
112 /* Set all SIMD register elements to 0 */
113 static gmx_inline gmx_simd_ref_pr
114 gmx_simd_ref_setzero_pr()
115 {
116     gmx_simd_ref_pr a;
117     int             i;
118
119     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
120     {
121         a.r[i] = 0.0;
122     }
123
124     return a;
125 }
126
127 static gmx_inline void
128 gmx_simd_ref_store_pr(real *dest, gmx_simd_ref_pr src)
129 {
130     int i;
131
132     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
133     {
134         dest[i] = src.r[i];
135     }
136 }
137
138 static gmx_inline gmx_simd_ref_pr
139 gmx_simd_ref_add_pr(gmx_simd_ref_pr a, gmx_simd_ref_pr b)
140 {
141     gmx_simd_ref_pr c;
142     int             i;
143
144     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
145     {
146         c.r[i] = a.r[i] + b.r[i];
147     }
148
149     return c;
150 }
151
152 static gmx_inline gmx_simd_ref_pr
153 gmx_simd_ref_sub_pr(gmx_simd_ref_pr a, gmx_simd_ref_pr b)
154 {
155     gmx_simd_ref_pr c;
156     int             i;
157
158     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
159     {
160         c.r[i] = a.r[i] - b.r[i];
161     }
162
163     return c;
164 }
165
166 static gmx_inline gmx_simd_ref_pr
167 gmx_simd_ref_mul_pr(gmx_simd_ref_pr a, gmx_simd_ref_pr b)
168 {
169     gmx_simd_ref_pr c;
170     int             i;
171
172     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
173     {
174         c.r[i] = a.r[i]*b.r[i];
175     }
176
177     return c;
178 }
179
180 static gmx_inline gmx_simd_ref_pr
181 gmx_simd_ref_madd_pr(gmx_simd_ref_pr a, gmx_simd_ref_pr b, gmx_simd_ref_pr c)
182 {
183     gmx_simd_ref_pr d;
184     int             i;
185
186     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
187     {
188         d.r[i] = a.r[i]*b.r[i] + c.r[i];
189     }
190
191     return d;
192 }
193
194 static gmx_inline gmx_simd_ref_pr
195 gmx_simd_ref_nmsub_pr(gmx_simd_ref_pr a, gmx_simd_ref_pr b, gmx_simd_ref_pr c)
196 {
197     gmx_simd_ref_pr d;
198     int             i;
199
200     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
201     {
202         d.r[i] = -a.r[i]*b.r[i] + c.r[i];
203     }
204
205     return d;
206 }
207
208 static gmx_inline gmx_simd_ref_pr
209 gmx_simd_ref_max_pr(gmx_simd_ref_pr a, gmx_simd_ref_pr b)
210 {
211     gmx_simd_ref_pr c;
212     int             i;
213
214     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
215     {
216         c.r[i] = (a.r[i] >= b.r[i] ? a.r[i] : b.r[i]);
217     }
218
219     return c;
220 }
221
222 static gmx_inline gmx_simd_ref_pr
223 gmx_simd_ref_blendzero_pr(gmx_simd_ref_pr a, gmx_simd_ref_pb b)
224 {
225     gmx_simd_ref_pr c;
226     int             i;
227
228     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
229     {
230         c.r[i] = (b.r[i] ? a.r[i] : 0.0);
231     }
232
233     return c;
234 }
235
236 /* Note that this reference implementation rounds away from zero,
237  * whereas most SIMD intrinsics will round to nearest even. Since this
238  * function is only used for periodic image calculations, the rounding
239  * of mantissas close to 0.5 is irrelevant, except in testing. This
240  * could be fixed by using rint/rintf, but the bigger problem is that
241  * MSVC does not support full C99, and none of the round or rint
242  * functions are defined. It's much easier to approximately implement
243  * round() than rint(), so we do that and hope we never get bitten in
244  * testing. (Thanks, Microsoft.)
245  */
246 static gmx_inline gmx_simd_ref_pr
247 gmx_simd_ref_round_pr(gmx_simd_ref_pr a)
248 {
249     gmx_simd_ref_pr b;
250     int             i;
251
252     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
253     {
254 #ifdef _MSC_VER
255         int temp = (a.r[i] >= 0.)
256             ? (a.r[i] + 0.5)
257             : (a.r[i] - 0.5);
258         b.r[i] = (real) temp;
259 #elif defined GMX_DOUBLE
260         b.r[i] = round(a.r[i]);
261 #else
262         b.r[i] = roundf(a.r[i]);
263 #endif
264     }
265
266     return b;
267 }
268
269 /* Not required, only used to speed up the nbnxn tabulated PME kernels */
270 static gmx_inline gmx_simd_ref_pr
271 gmx_simd_ref_floor_pr(gmx_simd_ref_pr a)
272 {
273     gmx_simd_ref_pr b;
274     int             i;
275
276     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
277     {
278 #ifdef GMX_DOUBLE
279         b.r[i] = floor(a.r[i]);
280 #else
281         b.r[i] = floorf(a.r[i]);
282 #endif
283     }
284
285     return b;
286 }
287
288 /* Not required, only used when blendv is faster than comparison */
289 static gmx_inline gmx_simd_ref_pr
290 gmx_simd_ref_blendv_pr(gmx_simd_ref_pr a, gmx_simd_ref_pr b, gmx_simd_ref_pr c)
291 {
292     gmx_simd_ref_pr d;
293     int             i;
294
295     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
296     {
297         d.r[i] = (c.r[i] >= 0) ? a.r[i] : b.r[i];
298     }
299
300     return d;
301 }
302
303 /* Copy the sign of a to b, assumes b >= 0 for efficiency */
304 static gmx_inline gmx_simd_ref_pr
305 gmx_simd_ref_cpsgn_nonneg_pr(gmx_simd_ref_pr a, gmx_simd_ref_pr b)
306 {
307     gmx_simd_ref_pr c;
308     int             i;
309
310     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
311     {
312         c.r[i] = (a.r[i] >= 0) ? b.r[i] : -b.r[i];
313     }
314
315     return c;
316 }
317
318 /* Very specific operation required in the non-bonded kernels */
319 static gmx_inline gmx_simd_ref_pr
320 gmx_simd_ref_masknot_add_pr(gmx_simd_ref_pb a, gmx_simd_ref_pr b, gmx_simd_ref_pr c)
321 {
322     gmx_simd_ref_pr d;
323     int             i;
324
325     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
326     {
327         d.r[i] = a.r[i] ? b.r[i] : b.r[i] + c.r[i];
328     }
329
330     return d;
331 }
332
333 /* Comparison */
334 static gmx_inline gmx_simd_ref_pb
335 gmx_simd_ref_cmplt_pr(gmx_simd_ref_pr a, gmx_simd_ref_pr b)
336 {
337     gmx_simd_ref_pb c;
338     int             i;
339
340     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
341     {
342         c.r[i] = (a.r[i] < b.r[i]);
343     }
344
345     return c;
346 }
347
348 /* Logical AND on SIMD booleans. */
349 static gmx_inline gmx_simd_ref_pb
350 gmx_simd_ref_and_pb(gmx_simd_ref_pb a, gmx_simd_ref_pb b)
351 {
352     gmx_simd_ref_pb c;
353     int             i;
354
355     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
356     {
357         c.r[i] = (a.r[i] && b.r[i]);
358     }
359
360     return c;
361 }
362
363 /* Logical OR on SIMD booleans. */
364 static gmx_inline gmx_simd_ref_pb
365 gmx_simd_ref_or_pb(gmx_simd_ref_pb a, gmx_simd_ref_pb b)
366 {
367     gmx_simd_ref_pb c;
368     int             i;
369
370     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
371     {
372         c.r[i] = (a.r[i] || b.r[i]);
373     }
374
375     return c;
376 }
377
378 /* Returns a single int (0/1) which tells if any of the booleans is True */
379 static gmx_inline int
380 gmx_simd_ref_anytrue_pb(gmx_simd_ref_pb a)
381 {
382     int anytrue;
383     int i;
384
385     anytrue = 0;
386     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
387     {
388         if (a.r[i])
389         {
390             anytrue = 1;
391         }
392     }
393
394     return anytrue;
395 }
396
397 /* Conversions only used for PME table lookup */
398 static gmx_inline gmx_simd_ref_epi32
399 gmx_simd_ref_cvttpr_epi32(gmx_simd_ref_pr a)
400 {
401     gmx_simd_ref_epi32 b;
402     int                i;
403
404     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
405     {
406         b.r[i] = (int)a.r[i];
407     }
408
409     return b;
410 };
411
412 /* These two function only need to be approximate, Newton-Raphson iteration
413  * is used for full accuracy in gmx_simd_invsqrt_r and gmx_simd_inv_r.
414  */
415 static gmx_inline gmx_simd_ref_pr
416 gmx_simd_ref_rsqrt_pr(gmx_simd_ref_pr a)
417 {
418     gmx_simd_ref_pr b;
419     int             i;
420
421     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
422     {
423 #ifdef GMX_DOUBLE
424         b.r[i] = 1.0/sqrt(a.r[i]);
425 #else
426         b.r[i] = 1.0/sqrtf(a.r[i]);
427 #endif
428     }
429
430     return b;
431 };
432
433 static gmx_inline gmx_simd_ref_pr
434 gmx_simd_ref_rcp_pr(gmx_simd_ref_pr a)
435 {
436     gmx_simd_ref_pr b;
437     int             i;
438
439     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
440     {
441         b.r[i] = 1.0/a.r[i];
442     }
443
444     return b;
445 };
446
447 static gmx_inline gmx_simd_ref_pr
448 gmx_simd_ref_exp_pr(gmx_simd_ref_pr a)
449 {
450     gmx_simd_ref_pr b;
451     int             i;
452
453     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
454     {
455 #ifdef GMX_DOUBLE
456         b.r[i] = exp(a.r[i]);
457 #else
458         b.r[i] = expf(a.r[i]);
459 #endif
460     }
461
462     return b;
463 };
464
465 static gmx_inline gmx_simd_ref_pr
466 gmx_simd_ref_sqrt_pr(gmx_simd_ref_pr a)
467 {
468     gmx_simd_ref_pr b;
469     int             i;
470
471     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
472     {
473 #ifdef GMX_DOUBLE
474         b.r[i] = sqrt(a.r[i]);
475 #else
476         b.r[i] = sqrtf(a.r[i]);
477 #endif
478     }
479
480     return b;
481 }
482
483 static gmx_inline int
484 gmx_simd_ref_sincos_pr(gmx_simd_ref_pr a,
485                        gmx_simd_ref_pr *s, gmx_simd_ref_pr *c)
486 {
487     int i;
488
489     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
490     {
491         s->r[i] = sin(a.r[i]);
492         c->r[i] = cos(a.r[i]);
493     }
494
495     return 0;
496 }
497
498 static gmx_inline gmx_simd_ref_pr
499 gmx_simd_ref_acos_pr(gmx_simd_ref_pr a)
500 {
501     gmx_simd_ref_pr b;
502     int             i;
503
504     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
505     {
506         b.r[i] = acos(a.r[i]);
507     }
508
509     return b;
510 }
511
512 static gmx_inline gmx_simd_ref_pr
513 gmx_simd_ref_atan2_pr(gmx_simd_ref_pr a, gmx_simd_ref_pr b)
514 {
515     gmx_simd_ref_pr c;
516     int             i;
517
518     for (i = 0; i < GMX_SIMD_REF_WIDTH; i++)
519     {
520         c.r[i] = atan2(a.r[i], b.r[i]);
521     }
522
523     return c;
524 }
525
526 #endif