Remove support for Intel classic compiler
[alexxy/gromacs.git] / src / gromacs / simd / tests / simd_floatingpoint_util.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2015,2017,2018,2019,2020,2021, 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 #include "gmxpre.h"
36
37 #include <numeric>
38
39 #include "gromacs/simd/simd.h"
40 #include "gromacs/utility/alignedallocator.h"
41 #include "gromacs/utility/basedefinitions.h"
42
43 #include "testutils/testasserts.h"
44
45 #include "simd.h"
46
47 namespace gmx
48 {
49 namespace test
50 {
51 namespace
52 {
53
54 /*! \cond internal */
55 /*! \addtogroup module_simd */
56 /*! \{ */
57
58 #if GMX_SIMD_HAVE_REAL
59
60 /*! \brief Test fixture for higher-level floating-point utility functions.
61  *
62  * Inherit from main SimdTest, add code to generate aligned memory and data.
63  */
64 class SimdFloatingpointUtilTest : public SimdTest
65 {
66 public:
67     SimdFloatingpointUtilTest()
68     {
69         // Resize vectors to get the amount of memory we need
70         integerMemory_.resize(GMX_SIMD_REAL_WIDTH);
71
72         // The total memory we allocate corresponds to two work arrays
73         // and 4 values each of GMX_SIMD_REAL_WIDTH.
74         realMemory_.resize(2 * s_workMemSize_ + 4 * GMX_SIMD_REAL_WIDTH);
75
76         offset_ = integerMemory_.data();
77         val0_   = realMemory_.data();
78         val1_   = val0_ + GMX_SIMD_REAL_WIDTH;
79         val2_   = val1_ + GMX_SIMD_REAL_WIDTH;
80         val3_   = val2_ + GMX_SIMD_REAL_WIDTH;
81         mem0_   = val3_ + GMX_SIMD_REAL_WIDTH;
82         mem1_   = mem0_ + s_workMemSize_;
83
84         // Set default values for offset and variables val0_ through val3_
85         // We cannot fill mem_ here since those values depend on the test.
86         for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
87         {
88             // Use every third point to avoid a continguous access pattern
89             offset_[i] = 3 * i;
90             // Multiply numbers by 1+100*GMX_REAL_EPS ensures some low bits are
91             // set too, so the tests make sure we read all bits correctly.
92             val0_[i] = (i) * (1.0 + 100 * GMX_REAL_EPS);
93             val1_[i] = (i + 0.1) * (1.0 + 100 * GMX_REAL_EPS);
94             val2_[i] = (i + 0.2) * (1.0 + 100 * GMX_REAL_EPS);
95             val3_[i] = (i + 0.3) * (1.0 + 100 * GMX_REAL_EPS);
96         }
97     }
98
99 protected:
100     //! \brief Size of memory work buffers
101     //
102     // To have a somewhat odd access pattern, we use every
103     // third entry, so the largest value of offset_[i] is 3*GMX_SIMD_REAL_WIDTH.
104     // Then we also allow alignments up to 16, which means the largest index in mem0_[]
105     // that we might access is 16*3*GMX_SIMD_REAL_WIDTH+3.
106     static const std::size_t s_workMemSize_ = 16 * 3 * GMX_SIMD_REAL_WIDTH + 4;
107
108     std::vector<int, AlignedAllocator<int>>   integerMemory_; //!< Aligned integer memory
109     std::vector<real, AlignedAllocator<real>> realMemory_;    //!< Aligned real memory
110
111     int*  offset_; //!< Pointer to offset indices, aligned memory
112     real* val0_;   //!< Pointer to GMX_SIMD_REAL_WIDTH values, aligned
113     real* val1_;   //!< Pointer to GMX_SIMD_REAL_WIDTH values, aligned
114     real* val2_;   //!< Pointer to GMX_SIMD_REAL_WIDTH values, aligned
115     real* val3_;   //!< Pointer to GMX_SIMD_REAL_WIDTH values, aligned
116
117     real* mem0_; //!< Pointer to aligned memory, s_workMemSize real values
118     real* mem1_; //!< Pointer to aligned memory, s_workMemSize real values
119 };
120
121
122 TEST_F(SimdFloatingpointUtilTest, gatherLoadTranspose4)
123 {
124     SimdReal  v0, v1, v2, v3;
125     SimdReal  ref0, ref1, ref2, ref3;
126     const int nalign                = 3;
127     int       alignmentList[nalign] = { 4, 8, 12 };
128     int       i, j, align;
129
130     for (i = 0; i < nalign; i++)
131     {
132         align = alignmentList[i];
133         for (j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
134         {
135             mem0_[align * offset_[j]]     = val0_[j];
136             mem0_[align * offset_[j] + 1] = val1_[j];
137             mem0_[align * offset_[j] + 2] = val2_[j];
138             mem0_[align * offset_[j] + 3] = val3_[j];
139         }
140
141         ref0 = load<SimdReal>(val0_);
142         ref1 = load<SimdReal>(val1_);
143         ref2 = load<SimdReal>(val2_);
144         ref3 = load<SimdReal>(val3_);
145
146         if (align == 4)
147         {
148             gatherLoadTranspose<4>(mem0_, offset_, &v0, &v1, &v2, &v3);
149         }
150         else if (align == 8)
151         {
152             gatherLoadTranspose<8>(mem0_, offset_, &v0, &v1, &v2, &v3);
153         }
154         else if (align == 12)
155         {
156             gatherLoadTranspose<12>(mem0_, offset_, &v0, &v1, &v2, &v3);
157         }
158         else
159         {
160             FAIL();
161         }
162
163         GMX_EXPECT_SIMD_REAL_EQ(ref0, v0);
164         GMX_EXPECT_SIMD_REAL_EQ(ref1, v1);
165         GMX_EXPECT_SIMD_REAL_EQ(ref2, v2);
166         GMX_EXPECT_SIMD_REAL_EQ(ref3, v3);
167     }
168 }
169
170 TEST_F(SimdFloatingpointUtilTest, gatherLoadTranspose2)
171 {
172     SimdReal  v0, v1;
173     SimdReal  ref0, ref1;
174     const int nalign                = 3;
175     int       alignmentList[nalign] = { 2, 4, c_simdBestPairAlignment };
176     int       i, j, align;
177
178     EXPECT_TRUE(c_simdBestPairAlignment <= GMX_SIMD_REAL_WIDTH);
179
180     for (i = 0; i < nalign; i++)
181     {
182         align = alignmentList[i];
183         for (j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
184         {
185             mem0_[align * offset_[j]]     = val0_[j];
186             mem0_[align * offset_[j] + 1] = val1_[j];
187         }
188
189         ref0 = load<SimdReal>(val0_);
190         ref1 = load<SimdReal>(val1_);
191
192         if (align == 2)
193         {
194             gatherLoadTranspose<2>(mem0_, offset_, &v0, &v1);
195         }
196         else if (align == 4)
197         {
198             gatherLoadTranspose<4>(mem0_, offset_, &v0, &v1);
199         }
200         else if (align == c_simdBestPairAlignment)
201         {
202             gatherLoadTranspose<c_simdBestPairAlignment>(mem0_, offset_, &v0, &v1);
203         }
204         else
205         {
206             FAIL();
207         }
208
209         GMX_EXPECT_SIMD_REAL_EQ(ref0, v0);
210         GMX_EXPECT_SIMD_REAL_EQ(ref1, v1);
211     }
212 }
213
214 TEST_F(SimdFloatingpointUtilTest, gatherLoadUTranspose3)
215 {
216     SimdReal  v0, v1, v2;
217     SimdReal  ref0, ref1, ref2;
218     const int nalign                = 2;
219     int       alignmentList[nalign] = { 3, 4 };
220     int       i, j, align;
221
222     for (i = 0; i < nalign; i++)
223     {
224         align = alignmentList[i];
225         for (j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
226         {
227             mem0_[align * offset_[j]]     = val0_[j];
228             mem0_[align * offset_[j] + 1] = val1_[j];
229             mem0_[align * offset_[j] + 2] = val2_[j];
230         }
231
232         ref0 = load<SimdReal>(val0_);
233         ref1 = load<SimdReal>(val1_);
234         ref2 = load<SimdReal>(val2_);
235
236         if (align == 3)
237         {
238             gatherLoadUTranspose<3>(mem0_, offset_, &v0, &v1, &v2);
239         }
240         else if (align == 4)
241         {
242             gatherLoadUTranspose<4>(mem0_, offset_, &v0, &v1, &v2);
243         }
244         else
245         {
246             FAIL();
247         }
248
249         GMX_EXPECT_SIMD_REAL_EQ(ref0, v0);
250         GMX_EXPECT_SIMD_REAL_EQ(ref1, v1);
251         GMX_EXPECT_SIMD_REAL_EQ(ref2, v2);
252     }
253 }
254
255 TEST_F(SimdFloatingpointUtilTest, transposeScatterStoreU3)
256 {
257     SimdReal               v0, v1, v2;
258     real                   refmem[s_workMemSize_];
259     const int              nalign                = 2;
260     int                    alignmentList[nalign] = { 3, 4 };
261     int                    i, align;
262     FloatingPointTolerance tolerance(defaultRealTolerance());
263
264     for (i = 0; i < nalign; i++)
265     {
266         align = alignmentList[i];
267
268         // Set test and reference memory to background value
269         for (std::size_t j = 0; j < s_workMemSize_; j++)
270         {
271             // Multiply by 1+100*eps to make sure low bits are also used
272             mem0_[j] = refmem[j] = (1000.0 + j) * (1.0 + 100 * GMX_REAL_EPS);
273         }
274
275         for (std::size_t j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
276         {
277             // set values in _reference_ memory (we will then test with mem0_, and compare)
278             refmem[align * offset_[j]]     = val0_[j];
279             refmem[align * offset_[j] + 1] = val1_[j];
280             refmem[align * offset_[j] + 2] = val2_[j];
281         }
282
283         v0 = load<SimdReal>(val0_);
284         v1 = load<SimdReal>(val1_);
285         v2 = load<SimdReal>(val2_);
286
287         if (align == 3)
288         {
289             transposeScatterStoreU<3>(mem0_, offset_, v0, v1, v2);
290         }
291         else if (align == 4)
292         {
293             transposeScatterStoreU<4>(mem0_, offset_, v0, v1, v2);
294         }
295         else
296         {
297             FAIL();
298         }
299
300         for (std::size_t j = 0; j < s_workMemSize_; j++)
301         {
302             EXPECT_REAL_EQ_TOL(refmem[j], mem0_[j], tolerance);
303         }
304     }
305 }
306
307 TEST_F(SimdFloatingpointUtilTest, transposeScatterIncrU3)
308 {
309     SimdReal               v0, v1, v2;
310     real                   refmem[s_workMemSize_];
311     const int              nalign                = 2;
312     int                    alignmentList[nalign] = { 3, 4 };
313     int                    i, align;
314     FloatingPointTolerance tolerance(defaultRealTolerance());
315
316     for (i = 0; i < nalign; i++)
317     {
318         align = alignmentList[i];
319
320         // Set test and reference memory to background value
321         for (std::size_t j = 0; j < s_workMemSize_; j++)
322         {
323             // Multiply by 1+100*eps to make sure low bits are also used
324             mem0_[j] = refmem[j] = (1000.0 + j) * (1.0 + 100 * GMX_REAL_EPS);
325         }
326
327         for (std::size_t j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
328         {
329             // Add values to _reference_ memory (we will then test with mem0_, and compare)
330             refmem[align * offset_[j]] += val0_[j];
331             refmem[align * offset_[j] + 1] += val1_[j];
332             refmem[align * offset_[j] + 2] += val2_[j];
333         }
334
335         v0 = load<SimdReal>(val0_);
336         v1 = load<SimdReal>(val1_);
337         v2 = load<SimdReal>(val2_);
338
339         if (align == 3)
340         {
341             transposeScatterIncrU<3>(mem0_, offset_, v0, v1, v2);
342         }
343         else if (align == 4)
344         {
345             transposeScatterIncrU<4>(mem0_, offset_, v0, v1, v2);
346         }
347         else
348         {
349             FAIL();
350         }
351
352         for (std::size_t j = 0; j < s_workMemSize_; j++)
353         {
354             EXPECT_REAL_EQ_TOL(refmem[j], mem0_[j], tolerance);
355         }
356     }
357 }
358
359 TEST_F(SimdFloatingpointUtilTest, transposeScatterIncrU3Overlapping)
360 {
361     SimdReal               v0, v1, v2;
362     real                   refmem[s_workMemSize_];
363     FloatingPointTolerance tolerance(defaultRealTolerance());
364
365     // Alter offset_ to make all entries point to the same (first) value, so all entries will overlap
366     for (std::size_t j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
367     {
368         offset_[j] = 0;
369     }
370
371     // Set test and reference memory to background value
372     for (std::size_t j = 0; j < s_workMemSize_; j++)
373     {
374         // Multiply by 1+100*eps to make sure low bits are also used
375         mem0_[j] = refmem[j] = (1000.0 + j) * (1.0 + 100 * GMX_REAL_EPS);
376     }
377
378     for (std::size_t j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
379     {
380         // Add values to _reference_ memory (we will then test with mem0_, and compare)
381         refmem[3 * offset_[j]] += val0_[j];
382         refmem[3 * offset_[j] + 1] += val1_[j];
383         refmem[3 * offset_[j] + 2] += val2_[j];
384     }
385
386     v0 = load<SimdReal>(val0_);
387     v1 = load<SimdReal>(val1_);
388     v2 = load<SimdReal>(val2_);
389
390     transposeScatterIncrU<3>(mem0_, offset_, v0, v1, v2);
391
392     for (std::size_t j = 0; j < s_workMemSize_; j++)
393     {
394         EXPECT_REAL_EQ_TOL(refmem[j], mem0_[j], tolerance);
395     }
396 }
397
398 TEST_F(SimdFloatingpointUtilTest, transposeScatterDecrU3)
399 {
400     SimdReal               v0, v1, v2;
401     real                   refmem[s_workMemSize_];
402     const int              nalign                = 2;
403     int                    alignmentList[nalign] = { 3, 4 };
404     int                    i, align;
405     FloatingPointTolerance tolerance(defaultRealTolerance());
406
407     for (i = 0; i < nalign; i++)
408     {
409         align = alignmentList[i];
410
411         // Set test and reference memory to background value
412         for (std::size_t j = 0; j < s_workMemSize_; j++)
413         {
414             // Multiply by 1+100*eps to make sure low bits are also used
415             mem0_[j] = refmem[j] = (1000.0 + j) * (1.0 + 100 * GMX_REAL_EPS);
416         }
417
418         for (std::size_t j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
419         {
420             // Subtract values from _reference_ memory (we will then test with mem0_, and compare)
421             refmem[align * offset_[j]] -= val0_[j];
422             refmem[align * offset_[j] + 1] -= val1_[j];
423             refmem[align * offset_[j] + 2] -= val2_[j];
424         }
425
426         v0 = load<SimdReal>(val0_);
427         v1 = load<SimdReal>(val1_);
428         v2 = load<SimdReal>(val2_);
429
430         if (align == 3)
431         {
432             transposeScatterDecrU<3>(mem0_, offset_, v0, v1, v2);
433         }
434         else if (align == 4)
435         {
436             transposeScatterDecrU<4>(mem0_, offset_, v0, v1, v2);
437         }
438         else
439         {
440             FAIL();
441         }
442
443         for (std::size_t j = 0; j < s_workMemSize_; j++)
444         {
445             EXPECT_REAL_EQ_TOL(refmem[j], mem0_[j], tolerance);
446         }
447     }
448 }
449
450 TEST_F(SimdFloatingpointUtilTest, transposeScatterDecrU3Overlapping)
451 {
452     SimdReal               v0, v1, v2;
453     real                   refmem[s_workMemSize_];
454     FloatingPointTolerance tolerance(defaultRealTolerance());
455
456     // Alter offset_ to make all entries point to the same (first) value, so all entries will overlap
457     for (std::size_t j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
458     {
459         offset_[j] = 0;
460     }
461
462     // Set test and reference memory to background value
463     for (std::size_t j = 0; j < s_workMemSize_; j++)
464     {
465         // Multiply by 1+100*eps to make sure low bits are also used
466         mem0_[j] = refmem[j] = (1000.0 + j) * (1.0 + 100 * GMX_REAL_EPS);
467     }
468
469     for (std::size_t j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
470     {
471         // Subtract values from _reference_ memory (we will then test with mem0_, and compare)
472         refmem[3 * offset_[j]] -= val0_[j];
473         refmem[3 * offset_[j] + 1] -= val1_[j];
474         refmem[3 * offset_[j] + 2] -= val2_[j];
475     }
476
477     v0 = load<SimdReal>(val0_);
478     v1 = load<SimdReal>(val1_);
479     v2 = load<SimdReal>(val2_);
480
481     transposeScatterDecrU<3>(mem0_, offset_, v0, v1, v2);
482
483     for (std::size_t j = 0; j < s_workMemSize_; j++)
484     {
485         EXPECT_REAL_EQ_TOL(refmem[j], mem0_[j], tolerance);
486     }
487 }
488
489 TEST_F(SimdFloatingpointUtilTest, expandScalarsToTriplets)
490 {
491     SimdReal vs, v0, v1, v2;
492     int      i;
493
494     for (i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
495     {
496         mem0_[i] = i;
497     }
498
499     vs = load<SimdReal>(mem0_);
500
501     expandScalarsToTriplets(vs, &v0, &v1, &v2);
502
503     store(val0_, v0);
504     store(val1_, v1);
505     store(val2_, v2);
506
507     for (i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
508     {
509         EXPECT_EQ(i / 3, val0_[i]);
510         EXPECT_EQ((i + GMX_SIMD_REAL_WIDTH) / 3, val1_[i]);
511         EXPECT_EQ((i + 2 * GMX_SIMD_REAL_WIDTH) / 3, val2_[i]);
512     }
513 }
514
515
516 TEST_F(SimdFloatingpointUtilTest, gatherLoadBySimdIntTranspose4)
517 {
518     SimdReal  v0, v1, v2, v3;
519     SimdReal  ref0, ref1, ref2, ref3;
520     SimdInt32 simdoffset;
521     const int nalign                = 3;
522     int       alignmentList[nalign] = { 4, 8, 12 };
523     int       i, j, align;
524
525     for (i = 0; i < nalign; i++)
526     {
527         align = alignmentList[i];
528         for (j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
529         {
530             mem0_[align * offset_[j]]     = val0_[j];
531             mem0_[align * offset_[j] + 1] = val1_[j];
532             mem0_[align * offset_[j] + 2] = val2_[j];
533             mem0_[align * offset_[j] + 3] = val3_[j];
534         }
535
536         simdoffset = load<SimdInt32>(offset_);
537         ref0       = load<SimdReal>(val0_);
538         ref1       = load<SimdReal>(val1_);
539         ref2       = load<SimdReal>(val2_);
540         ref3       = load<SimdReal>(val3_);
541
542         if (align == 4)
543         {
544             gatherLoadBySimdIntTranspose<4>(mem0_, simdoffset, &v0, &v1, &v2, &v3);
545         }
546         else if (align == 8)
547         {
548             gatherLoadBySimdIntTranspose<8>(mem0_, simdoffset, &v0, &v1, &v2, &v3);
549         }
550         else if (align == 12)
551         {
552             gatherLoadBySimdIntTranspose<12>(mem0_, simdoffset, &v0, &v1, &v2, &v3);
553         }
554         else
555         {
556             FAIL();
557         }
558
559         GMX_EXPECT_SIMD_REAL_EQ(ref0, v0);
560         GMX_EXPECT_SIMD_REAL_EQ(ref1, v1);
561         GMX_EXPECT_SIMD_REAL_EQ(ref2, v2);
562         GMX_EXPECT_SIMD_REAL_EQ(ref3, v3);
563     }
564 }
565
566
567 TEST_F(SimdFloatingpointUtilTest, gatherLoadBySimdIntTranspose2)
568 {
569     SimdReal  v0, v1;
570     SimdReal  ref0, ref1;
571     SimdInt32 simdoffset;
572     const int nalign                = 3;
573     int       alignmentList[nalign] = { 4, 8, 12 };
574     int       i, j, align;
575
576     for (i = 0; i < nalign; i++)
577     {
578         align = alignmentList[i];
579         for (j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
580         {
581             mem0_[align * offset_[j]]     = val0_[j];
582             mem0_[align * offset_[j] + 1] = val1_[j];
583         }
584
585         simdoffset = load<SimdInt32>(offset_);
586         ref0       = load<SimdReal>(val0_);
587         ref1       = load<SimdReal>(val1_);
588
589         if (align == 4)
590         {
591             gatherLoadBySimdIntTranspose<4>(mem0_, simdoffset, &v0, &v1);
592         }
593         else if (align == 8)
594         {
595             gatherLoadBySimdIntTranspose<8>(mem0_, simdoffset, &v0, &v1);
596         }
597         else if (align == 12)
598         {
599             gatherLoadBySimdIntTranspose<12>(mem0_, simdoffset, &v0, &v1);
600         }
601         else
602         {
603             FAIL();
604         }
605
606         GMX_EXPECT_SIMD_REAL_EQ(ref0, v0);
607         GMX_EXPECT_SIMD_REAL_EQ(ref1, v1);
608     }
609 }
610
611 #    if GMX_SIMD_HAVE_GATHER_LOADU_BYSIMDINT_TRANSPOSE_REAL
612 TEST_F(SimdFloatingpointUtilTest, gatherLoadUBySimdIntTranspose2)
613 {
614     SimdReal  v0, v1;
615     SimdReal  ref0, ref1;
616     SimdInt32 simdoffset;
617     const int nalign                = 3;
618     int       alignmentList[nalign] = { 1, 3, 5 };
619     int       i, j, align;
620
621     for (i = 0; i < nalign; i++)
622     {
623         align = alignmentList[i];
624         for (j = 0; j < GMX_SIMD_REAL_WIDTH; j++)
625         {
626             mem0_[align * offset_[j]]     = val0_[j];
627             mem0_[align * offset_[j] + 1] = val1_[j];
628         }
629
630         simdoffset = load<SimdInt32>(offset_);
631         ref0       = load<SimdReal>(val0_);
632         ref1       = load<SimdReal>(val1_);
633
634         if (align == 1)
635         {
636             gatherLoadUBySimdIntTranspose<1>(mem0_, simdoffset, &v0, &v1);
637         }
638         else if (align == 3)
639         {
640             gatherLoadUBySimdIntTranspose<3>(mem0_, simdoffset, &v0, &v1);
641         }
642         else if (align == 5)
643         {
644             gatherLoadUBySimdIntTranspose<5>(mem0_, simdoffset, &v0, &v1);
645         }
646         else
647         {
648             FAIL();
649         }
650
651         GMX_EXPECT_SIMD_REAL_EQ(ref0, v0);
652         GMX_EXPECT_SIMD_REAL_EQ(ref1, v1);
653     }
654 }
655 #    endif // GMX_SIMD_HAVE_GATHER_LOADU_BYSIMDINT_TRANSPOSE_REAL
656
657 TEST_F(SimdFloatingpointUtilTest, reduceIncr4Sum)
658 {
659     int                    i;
660     SimdReal               v0, v1, v2, v3;
661     real                   sum0, sum1, sum2, sum3, tstsum;
662     FloatingPointTolerance tolerance(defaultRealTolerance());
663
664     v0 = load<SimdReal>(val0_);
665     v1 = load<SimdReal>(val1_);
666     v2 = load<SimdReal>(val2_);
667     v3 = load<SimdReal>(val3_);
668
669     sum0 = sum1 = sum2 = sum3 = 0;
670     for (i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
671     {
672         sum0 += val0_[i];
673         sum1 += val1_[i];
674         sum2 += val2_[i];
675         sum3 += val3_[i];
676     }
677
678     // Just put some numbers in memory so we check the addition is correct
679     mem0_[0] = c0;
680     mem0_[1] = c1;
681     mem0_[2] = c2;
682     mem0_[3] = c3;
683
684     tstsum = reduceIncr4ReturnSum(mem0_, v0, v1, v2, v3);
685
686     EXPECT_REAL_EQ_TOL(c0 + sum0, mem0_[0], tolerance);
687     EXPECT_REAL_EQ_TOL(c1 + sum1, mem0_[1], tolerance);
688     EXPECT_REAL_EQ_TOL(c2 + sum2, mem0_[2], tolerance);
689     EXPECT_REAL_EQ_TOL(c3 + sum3, mem0_[3], tolerance);
690
691     EXPECT_REAL_EQ_TOL(sum0 + sum1 + sum2 + sum3, tstsum, tolerance);
692 }
693
694 #    if GMX_SIMD_HAVE_HSIMD_UTIL_REAL
695
696 TEST_F(SimdFloatingpointUtilTest, loadDualHsimd)
697 {
698     SimdReal v0, v1;
699
700     // Point p to the upper half of val0_
701     real* p = val0_ + GMX_SIMD_REAL_WIDTH / 2;
702
703     v0 = load<SimdReal>(val0_);
704     v1 = loadDualHsimd(val0_, p);
705
706     GMX_EXPECT_SIMD_REAL_EQ(v0, v1);
707 }
708
709 TEST_F(SimdFloatingpointUtilTest, loadDuplicateHsimd)
710 {
711     SimdReal v0, v1;
712     int      i;
713     // Point p to the upper half of val0_
714     real* p = val0_ + GMX_SIMD_REAL_WIDTH / 2;
715     // Copy data so upper half is identical to lower
716     for (i = 0; i < GMX_SIMD_REAL_WIDTH / 2; i++)
717     {
718         p[i] = val0_[i];
719     }
720
721     v0 = load<SimdReal>(val0_);
722     v1 = loadDuplicateHsimd(val0_);
723
724     GMX_EXPECT_SIMD_REAL_EQ(v0, v1);
725 }
726
727
728 TEST_F(SimdFloatingpointUtilTest, loadU1DualHsimd)
729 {
730     SimdReal v0, v1;
731     int      i;
732     real     data[2] = { 1, 2 };
733
734     // Point p to the upper half of val0_
735     real* p = val0_ + GMX_SIMD_REAL_WIDTH / 2;
736     // Set all low elements to data[0], an high to data[1]
737     for (i = 0; i < GMX_SIMD_REAL_WIDTH / 2; i++)
738     {
739         val0_[i] = data[0];
740         p[i]     = data[1];
741     }
742
743     v0 = load<SimdReal>(val0_);
744     v1 = loadU1DualHsimd(data);
745
746     GMX_EXPECT_SIMD_REAL_EQ(v0, v1);
747 }
748
749
750 TEST_F(SimdFloatingpointUtilTest, storeDualHsimd)
751 {
752     SimdReal v0;
753     int      i;
754
755     // Point p to the upper half of val0_
756     real* p = val0_ + GMX_SIMD_REAL_WIDTH / 2;
757
758     v0 = load<SimdReal>(val2_);
759     storeDualHsimd(val0_, p, v0);
760
761     for (i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
762     {
763         EXPECT_EQ(val2_[i], val0_[i]);
764     }
765 }
766
767 TEST_F(SimdFloatingpointUtilTest, incrDualHsimd)
768 {
769     real     reference[GMX_SIMD_REAL_WIDTH];
770     SimdReal v0;
771
772     // Create reference values
773     for (std::size_t i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
774     {
775         reference[i] = val0_[i] + val2_[i];
776     }
777
778     // Point p to the upper half of val0_
779     real* p = val0_ + GMX_SIMD_REAL_WIDTH / 2;
780
781     v0 = load<SimdReal>(val2_);
782     incrDualHsimd(val0_, p, v0);
783
784     for (std::size_t i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
785     {
786         EXPECT_EQ(reference[i], val0_[i]);
787     }
788 }
789
790 TEST_F(SimdFloatingpointUtilTest, incrDualHsimdOverlapping)
791 {
792     real     reference[GMX_SIMD_REAL_WIDTH / 2];
793     SimdReal v0;
794
795     // Create reference values
796     for (std::size_t i = 0; i < GMX_SIMD_REAL_WIDTH / 2; i++)
797     {
798         reference[i] = val0_[i] + val2_[i] + val2_[GMX_SIMD_REAL_WIDTH / 2 + i];
799     }
800
801     v0 = load<SimdReal>(val2_);
802     incrDualHsimd(val0_, val0_, v0);
803
804     for (std::size_t i = 0; i < GMX_SIMD_REAL_WIDTH / 2; i++)
805     {
806         EXPECT_EQ(reference[i], val0_[i]);
807     }
808 }
809
810 TEST_F(SimdFloatingpointUtilTest, decr3Hsimd)
811 {
812     SimdReal               v0, v1, v2;
813     real                   ref[3 * GMX_SIMD_REAL_WIDTH / 2];
814     int                    i, j;
815     FloatingPointTolerance tolerance(defaultRealTolerance());
816
817     // Point p to the upper half of val1_
818     real* p = val1_ + GMX_SIMD_REAL_WIDTH / 2;
819     for (i = 0; i < GMX_SIMD_REAL_WIDTH / 2; i++)
820     {
821         ref[i] = val0_[i] - (val1_[i] + p[i]);
822     }
823     p = val2_ + GMX_SIMD_REAL_WIDTH / 2;
824     for (j = 0; j < GMX_SIMD_REAL_WIDTH / 2; i++, j++)
825     {
826         ref[i] = val0_[i] - (val2_[j] + p[j]);
827     }
828     p = val3_ + GMX_SIMD_REAL_WIDTH / 2;
829     for (j = 0; j < GMX_SIMD_REAL_WIDTH / 2; i++, j++)
830     {
831         ref[i] = val0_[i] - (val3_[j] + p[j]);
832     }
833
834     v0 = load<SimdReal>(val1_);
835     v1 = load<SimdReal>(val2_);
836     v2 = load<SimdReal>(val3_);
837     decr3Hsimd(val0_, v0, v1, v2);
838
839     for (i = 0; i < 3 * GMX_SIMD_REAL_WIDTH / 2; i++)
840     {
841         EXPECT_REAL_EQ_TOL(ref[i], val0_[i], tolerance);
842     }
843 }
844
845
846 TEST_F(SimdFloatingpointUtilTest, gatherLoadTranspose2Hsimd)
847 {
848     SimdReal v0, v1;
849     SimdReal ref0, ref1;
850
851     const int nalign                = 3;
852     int       alignmentList[nalign] = { 2, 4, c_simdBestPairAlignment };
853     int       i, j, align;
854
855     for (i = 0; i < nalign; i++)
856     {
857         align = alignmentList[i];
858         for (j = 0; j < GMX_SIMD_REAL_WIDTH / 2; j++)
859         {
860             // Use mem0_ as base for lower half
861             mem0_[align * offset_[j]]     = val0_[j];
862             mem0_[align * offset_[j] + 1] = val1_[j];
863             // Use mem1_ as base for upper half
864             mem1_[align * offset_[j]]     = val0_[GMX_SIMD_REAL_WIDTH / 2 + j];
865             mem1_[align * offset_[j] + 1] = val1_[GMX_SIMD_REAL_WIDTH / 2 + j];
866         }
867
868         ref0 = load<SimdReal>(val0_);
869         ref1 = load<SimdReal>(val1_);
870
871         if (align == 2)
872         {
873             gatherLoadTransposeHsimd<2>(mem0_, mem1_, offset_, &v0, &v1);
874         }
875         else if (align == 4)
876         {
877             gatherLoadTransposeHsimd<4>(mem0_, mem1_, offset_, &v0, &v1);
878         }
879         else if (align == c_simdBestPairAlignment)
880         {
881             gatherLoadTransposeHsimd<c_simdBestPairAlignment>(mem0_, mem1_, offset_, &v0, &v1);
882         }
883         else
884         {
885             FAIL();
886         }
887
888         GMX_EXPECT_SIMD_REAL_EQ(ref0, v0);
889         GMX_EXPECT_SIMD_REAL_EQ(ref1, v1);
890     }
891 }
892
893
894 TEST_F(SimdFloatingpointUtilTest, reduceIncr4SumHsimd)
895 {
896     int                    i;
897     SimdReal               v0, v1;
898     real                   sum0, sum1, sum2, sum3, tstsum;
899     FloatingPointTolerance tolerance(defaultRealTolerance());
900
901     // Use the half-SIMD storage in memory val0_ and val1_.
902     v0 = load<SimdReal>(val0_);
903     v1 = load<SimdReal>(val1_);
904
905     sum0 = sum1 = sum2 = sum3 = 0;
906     for (i = 0; i < GMX_SIMD_REAL_WIDTH / 2; i++)
907     {
908         sum0 += val0_[i];
909         sum1 += val0_[GMX_SIMD_REAL_WIDTH / 2 + i];
910         sum2 += val1_[i];
911         sum3 += val1_[GMX_SIMD_REAL_WIDTH / 2 + i];
912     }
913
914     // Just put some numbers in memory so we check the addition is correct
915     mem0_[0] = c0;
916     mem0_[1] = c1;
917     mem0_[2] = c2;
918     mem0_[3] = c3;
919
920     tstsum = reduceIncr4ReturnSumHsimd(mem0_, v0, v1);
921
922     EXPECT_REAL_EQ_TOL(c0 + sum0, mem0_[0], tolerance);
923     EXPECT_REAL_EQ_TOL(c1 + sum1, mem0_[1], tolerance);
924     EXPECT_REAL_EQ_TOL(c2 + sum2, mem0_[2], tolerance);
925     EXPECT_REAL_EQ_TOL(c3 + sum3, mem0_[3], tolerance);
926
927     EXPECT_REAL_EQ_TOL(sum0 + sum1 + sum2 + sum3, tstsum, tolerance);
928 }
929
930 #    endif // GMX_SIMD_HAVE_HSIMD_UTIL_REAL
931
932 // Test Currently doesn't work for GMX_SIMD_REAL_WIDTH<4. Should be fixed by having GMX_EXPECT_SIMD_REAL_EQ which works for both Simd and Simd4
933 #    if GMX_SIMD_HAVE_4NSIMD_UTIL_REAL && GMX_SIMD_REAL_WIDTH >= 4
934
935 TEST_F(SimdFloatingpointUtilTest, loadUNDuplicate4)
936 {
937     Simd4NReal v0, v1;
938     int        i;
939     real       data[GMX_SIMD_REAL_WIDTH / 4];
940     std::iota(data, data + GMX_SIMD_REAL_WIDTH / 4, 1);
941
942     for (i = 0; i < GMX_SIMD_REAL_WIDTH / 4; i++)
943     {
944         val0_[i * 4] = val0_[i * 4 + 1] = val0_[i * 4 + 2] = val0_[i * 4 + 3] = data[i];
945     }
946
947     v0 = load<Simd4NReal>(val0_);
948     v1 = loadUNDuplicate4(data);
949
950     GMX_EXPECT_SIMD_REAL_EQ(v0, v1);
951 }
952
953 TEST_F(SimdFloatingpointUtilTest, load4DuplicateN)
954 {
955     Simd4NReal v0, v1;
956     int        i;
957     real       data[4] = { 1, 2, 3, 4 };
958
959     for (i = 0; i < GMX_SIMD_REAL_WIDTH / 4; i++)
960     {
961         val0_[i * 4]     = data[0];
962         val0_[i * 4 + 1] = data[1];
963         val0_[i * 4 + 2] = data[2];
964         val0_[i * 4 + 3] = data[3];
965     }
966
967     v0 = load<Simd4NReal>(val0_);
968     v1 = load4DuplicateN(val0_);
969
970     GMX_EXPECT_SIMD_REAL_EQ(v0, v1);
971 }
972
973 TEST_F(SimdFloatingpointUtilTest, loadU4NOffset)
974 {
975     constexpr int offset  = 6; // non power of 2
976     constexpr int dataLen = 4 + offset * (GMX_SIMD_REAL_WIDTH / 4 - 1);
977     real          data[dataLen];
978     std::iota(data, data + dataLen, 1);
979
980     for (int i = 0; i < GMX_SIMD_REAL_WIDTH / 4; i++)
981     {
982         val0_[i * 4]     = data[0 + offset * i];
983         val0_[i * 4 + 1] = data[1 + offset * i];
984         val0_[i * 4 + 2] = data[2 + offset * i];
985         val0_[i * 4 + 3] = data[3 + offset * i];
986     }
987
988     const Simd4NReal v0 = load<Simd4NReal>(val0_);
989     const Simd4NReal v1 = loadU4NOffset(data, offset);
990
991     GMX_EXPECT_SIMD_REAL_EQ(v0, v1);
992 }
993
994 #    endif // GMX_SIMD_HAVE_4NSIMD_UTIL_REAL
995
996 #endif // GMX_SIMD_HAVE_REAL
997
998 /*! \} */
999 /*! \endcond */
1000
1001 } // namespace
1002 } // namespace test
1003 } // namespace gmx