Fix compilation issues with ARM SIMD
[alexxy/gromacs.git] / src / gromacs / simd / impl_arm_neon_asimd / impl_arm_neon_asimd_simd_double.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,2017, 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_IMPL_ARM_NEON_ASIMD_SIMD_DOUBLE_H
37 #define GMX_SIMD_IMPL_ARM_NEON_ASIMD_SIMD_DOUBLE_H
38
39 #include "config.h"
40
41 #include <cassert>
42
43 #include <arm_neon.h>
44
45 #include "impl_arm_neon_asimd_simd_float.h"
46
47 namespace gmx
48 {
49
50 class SimdDouble
51 {
52     public:
53         SimdDouble() {}
54
55         SimdDouble(double d) : simdInternal_(vdupq_n_f64(d)) {}
56
57         // Internal utility constructor to simplify return statements
58         SimdDouble(float64x2_t simd) : simdInternal_(simd) {}
59
60         float64x2_t  simdInternal_;
61 };
62
63 class SimdDInt32
64 {
65     public:
66         SimdDInt32() {}
67
68         SimdDInt32(std::int32_t i) : simdInternal_(vdup_n_s32(i)) {}
69
70         // Internal utility constructor to simplify return statements
71         SimdDInt32(int32x2_t simd) : simdInternal_(simd) {}
72
73         int32x2_t  simdInternal_;
74 };
75
76 class SimdDBool
77 {
78     public:
79         SimdDBool() {}
80
81         SimdDBool(bool b) : simdInternal_(vdupq_n_u64( b ? 0xFFFFFFFFFFFFFFFF : 0)) {}
82
83         // Internal utility constructor to simplify return statements
84         SimdDBool(uint64x2_t simd) : simdInternal_(simd) {}
85
86         uint64x2_t  simdInternal_;
87 };
88
89 class SimdDIBool
90 {
91     public:
92         SimdDIBool() {}
93
94         SimdDIBool(bool b) : simdInternal_(vdup_n_u32( b ? 0xFFFFFFFF : 0)) {}
95
96         // Internal utility constructor to simplify return statements
97         SimdDIBool(uint32x2_t simd) : simdInternal_(simd) {}
98
99         uint32x2_t  simdInternal_;
100 };
101
102 static inline SimdDouble gmx_simdcall
103 simdLoad(const double *m)
104 {
105     assert(std::size_t(m) % 16 == 0);
106     return {
107                vld1q_f64(m)
108     };
109 }
110
111 static inline void gmx_simdcall
112 store(double *m, SimdDouble a)
113 {
114     assert(std::size_t(m) % 16 == 0);
115     vst1q_f64(m, a.simdInternal_);
116 }
117
118 static inline SimdDouble gmx_simdcall
119 simdLoadU(const double *m)
120 {
121     return {
122                vld1q_f64(m)
123     };
124 }
125
126 static inline void gmx_simdcall
127 storeU(double *m, SimdDouble a)
128 {
129     vst1q_f64(m, a.simdInternal_);
130 }
131
132 static inline SimdDouble gmx_simdcall
133 setZeroD()
134 {
135     return {
136                vdupq_n_f64(0.0)
137     };
138 }
139
140 static inline SimdDInt32 gmx_simdcall
141 simdLoadDI(const std::int32_t * m)
142 {
143     assert(std::size_t(m) % 8 == 0);
144     return {
145                vld1_s32(m)
146     };
147 }
148
149 static inline void gmx_simdcall
150 store(std::int32_t * m, SimdDInt32 a)
151 {
152     assert(std::size_t(m) % 8 == 0);
153     vst1_s32(m, a.simdInternal_);
154 }
155
156 static inline SimdDInt32 gmx_simdcall
157 simdLoadUDI(const std::int32_t *m)
158 {
159     return {
160                vld1_s32(m)
161     };
162 }
163
164 static inline void gmx_simdcall
165 storeU(std::int32_t * m, SimdDInt32 a)
166 {
167     vst1_s32(m, a.simdInternal_);
168 }
169
170 static inline SimdDInt32 gmx_simdcall
171 setZeroDI()
172 {
173     return {
174                vdup_n_s32(0)
175     };
176 }
177
178 template<int index> gmx_simdcall
179 static inline std::int32_t
180 extract(SimdDInt32 a)
181 {
182     return vget_lane_s32(a.simdInternal_, index);
183 }
184
185 static inline SimdDouble gmx_simdcall
186 operator&(SimdDouble a, SimdDouble b)
187 {
188     return {
189                float64x2_t(vandq_s64(int64x2_t(a.simdInternal_), int64x2_t(b.simdInternal_)))
190     };
191 }
192
193 static inline SimdDouble gmx_simdcall
194 andNot(SimdDouble a, SimdDouble b)
195 {
196     return {
197                float64x2_t(vbicq_s64(int64x2_t(b.simdInternal_), int64x2_t(a.simdInternal_)))
198     };
199 }
200
201 static inline SimdDouble gmx_simdcall
202 operator|(SimdDouble a, SimdDouble b)
203 {
204     return {
205                float64x2_t(vorrq_s64(int64x2_t(a.simdInternal_), int64x2_t(b.simdInternal_)))
206     };
207 }
208
209 static inline SimdDouble gmx_simdcall
210 operator^(SimdDouble a, SimdDouble b)
211 {
212     return {
213                float64x2_t(veorq_s64(int64x2_t(a.simdInternal_), int64x2_t(b.simdInternal_)))
214     };
215 }
216
217 static inline SimdDouble gmx_simdcall
218 operator+(SimdDouble a, SimdDouble b)
219 {
220     return {
221                vaddq_f64(a.simdInternal_, b.simdInternal_)
222     };
223 }
224
225 static inline SimdDouble gmx_simdcall
226 operator-(SimdDouble a, SimdDouble b)
227 {
228     return {
229                vsubq_f64(a.simdInternal_, b.simdInternal_)
230     };
231 }
232
233 static inline SimdDouble gmx_simdcall
234 operator-(SimdDouble x)
235 {
236     return {
237                vnegq_f64(x.simdInternal_)
238     };
239 }
240
241 static inline SimdDouble gmx_simdcall
242 operator*(SimdDouble a, SimdDouble b)
243 {
244     return {
245                vmulq_f64(a.simdInternal_, b.simdInternal_)
246     };
247 }
248
249 static inline SimdDouble gmx_simdcall
250 fma(SimdDouble a, SimdDouble b, SimdDouble c)
251 {
252     return {
253                vfmaq_f64(c.simdInternal_, b.simdInternal_, a.simdInternal_)
254     };
255 }
256
257 static inline SimdDouble gmx_simdcall
258 fms(SimdDouble a, SimdDouble b, SimdDouble c)
259 {
260     return {
261                vnegq_f64(vfmsq_f64(c.simdInternal_, b.simdInternal_, a.simdInternal_))
262     };
263 }
264
265 static inline SimdDouble gmx_simdcall
266 fnma(SimdDouble a, SimdDouble b, SimdDouble c)
267 {
268     return {
269                vfmsq_f64(c.simdInternal_, b.simdInternal_, a.simdInternal_)
270     };
271 }
272
273 static inline SimdDouble gmx_simdcall
274 fnms(SimdDouble a, SimdDouble b, SimdDouble c)
275 {
276     return {
277                vnegq_f64(vfmaq_f64(c.simdInternal_, b.simdInternal_, a.simdInternal_))
278     };
279 }
280
281 static inline SimdDouble gmx_simdcall
282 rsqrt(SimdDouble x)
283 {
284     return {
285                vrsqrteq_f64(x.simdInternal_)
286     };
287 }
288
289 static inline SimdDouble gmx_simdcall
290 rsqrtIter(SimdDouble lu, SimdDouble x)
291 {
292     return {
293                vmulq_f64(lu.simdInternal_, vrsqrtsq_f64(vmulq_f64(lu.simdInternal_, lu.simdInternal_), x.simdInternal_))
294     };
295 }
296
297 static inline SimdDouble gmx_simdcall
298 rcp(SimdDouble x)
299 {
300     return {
301                vrecpeq_f64(x.simdInternal_)
302     };
303 }
304
305 static inline SimdDouble gmx_simdcall
306 rcpIter(SimdDouble lu, SimdDouble x)
307 {
308     return {
309                vmulq_f64(lu.simdInternal_, vrecpsq_f64(lu.simdInternal_, x.simdInternal_))
310     };
311 }
312
313 static inline SimdDouble gmx_simdcall
314 maskAdd(SimdDouble a, SimdDouble b, SimdDBool m)
315 {
316     float64x2_t addend = float64x2_t(vandq_u64(uint64x2_t(b.simdInternal_), m.simdInternal_));
317
318     return {
319                vaddq_f64(a.simdInternal_, addend)
320     };
321 }
322
323 static inline SimdDouble gmx_simdcall
324 maskzMul(SimdDouble a, SimdDouble b, SimdDBool m)
325 {
326     float64x2_t prod = vmulq_f64(a.simdInternal_, b.simdInternal_);
327     return {
328                float64x2_t(vandq_u64(uint64x2_t(prod), m.simdInternal_))
329     };
330 }
331
332 static inline SimdDouble gmx_simdcall
333 maskzFma(SimdDouble a, SimdDouble b, SimdDouble c, SimdDBool m)
334 {
335     float64x2_t prod = vfmaq_f64(c.simdInternal_, b.simdInternal_, a.simdInternal_);
336
337     return {
338                float64x2_t(vandq_u64(uint64x2_t(prod), m.simdInternal_))
339     };
340 }
341
342 static inline SimdDouble gmx_simdcall
343 maskzRsqrt(SimdDouble x, SimdDBool m)
344 {
345     // The result will always be correct since we mask the result with m, but
346     // for debug builds we also want to make sure not to generate FP exceptions
347 #ifndef NDEBUG
348     x.simdInternal_ = vbslq_f64(m.simdInternal_, x.simdInternal_, vdupq_n_f64(1.0f));
349 #endif
350     return {
351                float64x2_t(vandq_u64(uint64x2_t(vrsqrteq_f64(x.simdInternal_)), m.simdInternal_))
352     };
353 }
354
355 static inline SimdDouble gmx_simdcall
356 maskzRcp(SimdDouble x, SimdDBool m)
357 {
358     // The result will always be correct since we mask the result with m, but
359     // for debug builds we also want to make sure not to generate FP exceptions
360 #ifndef NDEBUG
361     x.simdInternal_ = vbslq_f64(m.simdInternal_, x.simdInternal_, vdupq_n_f64(1.0f));
362 #endif
363     return {
364                float64x2_t(vandq_u64(uint64x2_t(vrecpeq_f64(x.simdInternal_)), m.simdInternal_))
365     };
366 }
367
368 static inline SimdDouble gmx_simdcall
369 abs(SimdDouble x)
370 {
371     return {
372                vabsq_f64( x.simdInternal_ )
373     };
374 }
375
376 static inline SimdDouble gmx_simdcall
377 max(SimdDouble a, SimdDouble b)
378 {
379     return {
380                vmaxq_f64(a.simdInternal_, b.simdInternal_)
381     };
382 }
383
384 static inline SimdDouble gmx_simdcall
385 min(SimdDouble a, SimdDouble b)
386 {
387     return {
388                vminq_f64(a.simdInternal_, b.simdInternal_)
389     };
390 }
391
392 static inline SimdDouble gmx_simdcall
393 round(SimdDouble x)
394 {
395     return {
396                vrndnq_f64(x.simdInternal_)
397     };
398 }
399
400 static inline SimdDouble gmx_simdcall
401 trunc(SimdDouble x)
402 {
403     return {
404                vrndq_f64( x.simdInternal_ )
405     };
406 }
407
408 static inline SimdDouble
409 frexp(SimdDouble value, SimdDInt32 * exponent)
410 {
411     const float64x2_t exponentMask = float64x2_t( vdupq_n_s64(0x7FF0000000000000LL) );
412     const float64x2_t mantissaMask = float64x2_t( vdupq_n_s64(0x800FFFFFFFFFFFFFLL) );
413
414     const int64x2_t   exponentBias = vdupq_n_s64(1022); // add 1 to make our definition identical to frexp()
415     const float64x2_t half         = vdupq_n_f64(0.5);
416     int64x2_t         iExponent;
417
418     iExponent               = vandq_s64( int64x2_t(value.simdInternal_), int64x2_t(exponentMask) );
419     iExponent               = vsubq_s64(vshrq_n_s64(iExponent, 52), exponentBias);
420     exponent->simdInternal_ = vmovn_s64(iExponent);
421
422     return {
423                float64x2_t(vorrq_s64(vandq_s64(int64x2_t(value.simdInternal_), int64x2_t(mantissaMask)), int64x2_t(half)))
424     };
425 }
426
427 static inline SimdDouble
428 ldexp(SimdDouble value, SimdDInt32 exponent)
429 {
430     const int64x2_t  exponentBias = vdupq_n_s64(1023);
431     int64x2_t        iExponent;
432
433     iExponent = vmovl_s32(exponent.simdInternal_);
434     iExponent = vshlq_n_s64(vaddq_s64(iExponent, exponentBias), 52);
435
436     return {
437                vmulq_f64(value.simdInternal_, float64x2_t(iExponent))
438     };
439 }
440
441 static inline double gmx_simdcall
442 reduce(SimdDouble a)
443 {
444     float64x2_t b = vpaddq_f64(a.simdInternal_, a.simdInternal_);
445     return vgetq_lane_f64(b, 0);
446 }
447
448 static inline SimdDBool gmx_simdcall
449 operator==(SimdDouble a, SimdDouble b)
450 {
451     return {
452                vceqq_f64(a.simdInternal_, b.simdInternal_)
453     };
454 }
455
456 static inline SimdDBool gmx_simdcall
457 operator!=(SimdDouble a, SimdDouble b)
458 {
459     return {
460                vreinterpretq_u64_u32(vmvnq_u32(vreinterpretq_u32_u64(vceqq_f64(a.simdInternal_, b.simdInternal_))))
461     };
462 }
463
464 static inline SimdDBool gmx_simdcall
465 operator<(SimdDouble a, SimdDouble b)
466 {
467     return {
468                vcltq_f64(a.simdInternal_, b.simdInternal_)
469     };
470 }
471
472 static inline SimdDBool gmx_simdcall
473 operator<=(SimdDouble a, SimdDouble b)
474 {
475     return {
476                vcleq_f64(a.simdInternal_, b.simdInternal_)
477     };
478 }
479
480 static inline SimdDBool gmx_simdcall
481 testBits(SimdDouble a)
482 {
483     return {
484                vtstq_s64( int64x2_t(a.simdInternal_), int64x2_t(a.simdInternal_) )
485     };
486 }
487
488 static inline SimdDBool gmx_simdcall
489 operator&&(SimdDBool a, SimdDBool b)
490 {
491     return {
492                vandq_u64(a.simdInternal_, b.simdInternal_)
493     };
494 }
495
496 static inline SimdDBool gmx_simdcall
497 operator||(SimdDBool a, SimdDBool b)
498 {
499     return {
500                vorrq_u64(a.simdInternal_, b.simdInternal_)
501     };
502 }
503
504 static inline bool gmx_simdcall
505 anyTrue(SimdDBool a)
506 {
507     return (vmaxvq_u32((uint32x4_t)(a.simdInternal_)) != 0);
508 }
509
510 static inline SimdDouble gmx_simdcall
511 selectByMask(SimdDouble a, SimdDBool m)
512 {
513     return {
514         float64x2_t(vandq_u64(uint64x2_t(a.simdInternal_), m.simdInternal_))
515     };
516 }
517
518 static inline SimdDouble gmx_simdcall
519 selectByNotMask(SimdDouble a, SimdDBool m)
520 {
521     return {
522         float64x2_t(vbicq_u64(uint64x2_t(a.simdInternal_), m.simdInternal_))
523     };
524 }
525
526 static inline SimdDouble gmx_simdcall
527 blend(SimdDouble a, SimdDouble b, SimdDBool sel)
528 {
529     return {
530         vbslq_f64(sel.simdInternal_, b.simdInternal_, a.simdInternal_)
531     };
532 }
533
534 static inline SimdDInt32 gmx_simdcall
535 operator<<(SimdDInt32 a, int n)
536 {
537     return {
538         vshl_s32(a.simdInternal_, vdup_n_s32(n >= 32 ? 32 : n))
539     };
540 }
541
542 static inline SimdDInt32 gmx_simdcall
543 operator>>(SimdDInt32 a, int n)
544 {
545     return {
546         vshl_s32(a.simdInternal_, vdup_n_s32(n >= 32 ? -32 : -n))
547     };
548 }
549
550 static inline SimdDInt32 gmx_simdcall
551 operator&(SimdDInt32 a, SimdDInt32 b)
552 {
553     return {
554         vand_s32(a.simdInternal_, b.simdInternal_)
555     };
556 }
557
558 static inline SimdDInt32 gmx_simdcall
559 andNot(SimdDInt32 a, SimdDInt32 b)
560 {
561     return {
562         vbic_s32(b.simdInternal_, a.simdInternal_)
563     };
564 }
565
566 static inline SimdDInt32 gmx_simdcall
567 operator|(SimdDInt32 a, SimdDInt32 b)
568 {
569     return {
570         vorr_s32(a.simdInternal_, b.simdInternal_)
571     };
572 }
573
574 static inline SimdDInt32 gmx_simdcall
575 operator^(SimdDInt32 a, SimdDInt32 b)
576 {
577     return {
578         veor_s32(a.simdInternal_, b.simdInternal_)
579     };
580 }
581
582 static inline SimdDInt32 gmx_simdcall
583 operator+(SimdDInt32 a, SimdDInt32 b)
584 {
585     return {
586         vadd_s32(a.simdInternal_, b.simdInternal_)
587     };
588 }
589
590 static inline SimdDInt32 gmx_simdcall
591 operator-(SimdDInt32 a, SimdDInt32 b)
592 {
593     return {
594         vsub_s32(a.simdInternal_, b.simdInternal_)
595     };
596 }
597
598 static inline SimdDInt32 gmx_simdcall
599 operator*(SimdDInt32 a, SimdDInt32 b)
600 {
601     return {
602         vmul_s32(a.simdInternal_, b.simdInternal_)
603     };
604 }
605
606 static inline SimdDIBool gmx_simdcall
607 operator==(SimdDInt32 a, SimdDInt32 b)
608 {
609     return {
610         vceq_s32(a.simdInternal_, b.simdInternal_)
611     };
612 }
613
614 static inline SimdDIBool gmx_simdcall
615 testBits(SimdDInt32 a)
616 {
617     return {
618         vtst_s32( a.simdInternal_, a.simdInternal_)
619     };
620 }
621
622 static inline SimdDIBool gmx_simdcall
623 operator<(SimdDInt32 a, SimdDInt32 b)
624 {
625     return {
626         vclt_s32(a.simdInternal_, b.simdInternal_)
627     };
628 }
629
630 static inline SimdDIBool gmx_simdcall
631 operator&&(SimdDIBool a, SimdDIBool b)
632 {
633     return {
634         vand_u32(a.simdInternal_, b.simdInternal_)
635     };
636 }
637
638 static inline SimdDIBool gmx_simdcall
639 operator||(SimdDIBool a, SimdDIBool b)
640 {
641     return {
642         vorr_u32(a.simdInternal_, b.simdInternal_)
643     };
644 }
645
646 static inline bool gmx_simdcall
647 anyTrue(SimdDIBool a)
648 {
649     return (vmaxv_u32(a.simdInternal_) != 0);
650 }
651
652 static inline SimdDInt32 gmx_simdcall
653 selectByMask(SimdDInt32 a, SimdDIBool m)
654 {
655     return {
656         vand_s32(a.simdInternal_, vreinterpret_s32_u32(m.simdInternal_))
657     };
658 }
659
660 static inline SimdDInt32 gmx_simdcall
661 selectByNotMask(SimdDInt32 a, SimdDIBool m)
662 {
663     return {
664         vbic_s32(a.simdInternal_, vreinterpret_s32_u32(m.simdInternal_))
665     };
666 }
667
668 static inline SimdDInt32 gmx_simdcall
669 blend(SimdDInt32 a, SimdDInt32 b, SimdDIBool sel)
670 {
671     return {
672         vbsl_s32(sel.simdInternal_, b.simdInternal_, a.simdInternal_)
673     };
674 }
675
676 static inline SimdDInt32 gmx_simdcall
677 cvtR2I(SimdDouble a)
678 {
679     return {
680         vmovn_s64(vcvtnq_s64_f64(a.simdInternal_))
681     };
682 }
683
684 static inline SimdDInt32 gmx_simdcall
685 cvttR2I(SimdDouble a)
686 {
687     return {
688         vmovn_s64(vcvtq_s64_f64(a.simdInternal_))
689     };
690 }
691
692 static inline SimdDouble gmx_simdcall
693 cvtI2R(SimdDInt32 a)
694 {
695     return {
696         vcvtq_f64_s64(vmovl_s32(a.simdInternal_))
697     };
698 }
699
700 static inline SimdDIBool gmx_simdcall
701 cvtB2IB(SimdDBool a)
702 {
703     return {
704         vqmovn_u64(a.simdInternal_)
705     };
706 }
707
708 static inline SimdDBool gmx_simdcall
709 cvtIB2B(SimdDIBool a)
710 {
711     return {
712         vorrq_u64(vmovl_u32(a.simdInternal_), vshlq_n_u64(vmovl_u32(a.simdInternal_), 32))
713     };
714 }
715
716 static inline void gmx_simdcall
717 cvtF2DD(SimdFloat f, SimdDouble *d0, SimdDouble *d1)
718 {
719     d0->simdInternal_ = vcvt_f64_f32(vget_low_f32(f.simdInternal_));
720     d1->simdInternal_ = vcvt_high_f64_f32(f.simdInternal_);
721 }
722
723 static inline SimdFloat gmx_simdcall
724 cvtDD2F(SimdDouble d0, SimdDouble d1)
725 {
726     return {
727         vcvt_high_f32_f64(vcvt_f32_f64(d0.simdInternal_), d1.simdInternal_)
728     };
729 }
730
731 }      // namespace gmx
732
733 #endif // GMX_SIMD_IMPL_ARM_NEON_ASIMD_SIMD_DOUBLE_H