Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / simd / tests / simd_vector_operations.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 #include "gmxpre.h"
36
37 #include <math.h>
38
39 #include "gromacs/simd/simd.h"
40 #include "gromacs/simd/vector_operations.h"
41
42 #include "simd.h"
43
44 namespace gmx
45 {
46 namespace test
47 {
48 namespace
49 {
50
51 /*! \cond internal */
52 /*! \addtogroup module_simd */
53 /*! \{ */
54
55 #ifdef GMX_SIMD_HAVE_REAL
56
57 /*! \internal \brief Test fixture for vector operations tests (identical to the generic \ref SimdTest) */
58 typedef SimdTest SimdVectorOperationsTest;
59
60 TEST_F(SimdVectorOperationsTest, gmxSimdCalcRsqR)
61 {
62     gmx_simd_real_t simdX  = setSimdRealFrom3R(1, 2, 3);
63     gmx_simd_real_t simdY  = setSimdRealFrom3R(3, 0, 5);
64     gmx_simd_real_t simdZ  = setSimdRealFrom3R(4, 1, 8);
65     gmx_simd_real_t simdR2 = setSimdRealFrom3R(26, 5, 98);
66
67     setUlpTol(2);
68     GMX_EXPECT_SIMD_REAL_NEAR(simdR2, gmx_simd_calc_rsq_r(simdX, simdY, simdZ));
69 }
70
71 TEST_F(SimdVectorOperationsTest, gmxSimdIprodR)
72 {
73     gmx_simd_real_t aX    = setSimdRealFrom3R(1, 2, 3);
74     gmx_simd_real_t aY    = setSimdRealFrom3R(3, 0, 5);
75     gmx_simd_real_t aZ    = setSimdRealFrom3R(4, 1, 8);
76     gmx_simd_real_t bX    = setSimdRealFrom3R(8, 3, 6);
77     gmx_simd_real_t bY    = setSimdRealFrom3R(2, 3, 1);
78     gmx_simd_real_t bZ    = setSimdRealFrom3R(5, 7, 9);
79     gmx_simd_real_t iprod = setSimdRealFrom3R(34, 13, 95);
80
81     setUlpTol(2);
82     GMX_EXPECT_SIMD_REAL_NEAR(iprod, gmx_simd_iprod_r(aX, aY, aZ, bX, bY, bZ));
83 }
84
85 TEST_F(SimdVectorOperationsTest, gmxSimdNorm2R)
86 {
87     gmx_simd_real_t simdX     = setSimdRealFrom3R(1, 2, 3);
88     gmx_simd_real_t simdY     = setSimdRealFrom3R(3, 0, 5);
89     gmx_simd_real_t simdZ     = setSimdRealFrom3R(4, 1, 8);
90     gmx_simd_real_t simdNorm2 = setSimdRealFrom3R(26, 5, 98);
91
92     setUlpTol(2);
93     GMX_EXPECT_SIMD_REAL_NEAR(simdNorm2, gmx_simd_norm2_r(simdX, simdY, simdZ));
94 }
95
96 TEST_F(SimdVectorOperationsTest, gmxSimdCprodR)
97 {
98     gmx_simd_real_t aX    = setSimdRealFrom3R(1, 2, 3);
99     gmx_simd_real_t aY    = setSimdRealFrom3R(3, 0, 5);
100     gmx_simd_real_t aZ    = setSimdRealFrom3R(4, 1, 8);
101     gmx_simd_real_t bX    = setSimdRealFrom3R(8, 3, 6);
102     gmx_simd_real_t bY    = setSimdRealFrom3R(2, 3, 1);
103     gmx_simd_real_t bZ    = setSimdRealFrom3R(5, 7, 9);
104     gmx_simd_real_t refcX = setSimdRealFrom3R(7, -3, 37);
105     gmx_simd_real_t refcY = setSimdRealFrom3R(27, -11, 21);
106     gmx_simd_real_t refcZ = setSimdRealFrom3R(-22, 6, -27);
107     gmx_simd_real_t cX, cY, cZ;
108     gmx_simd_cprod_r(aX, aY, aZ, bX, bY, bZ, &cX, &cY, &cZ);
109
110     setUlpTol(2);
111     GMX_EXPECT_SIMD_REAL_NEAR(refcX, cX);
112     GMX_EXPECT_SIMD_REAL_NEAR(refcY, cY);
113     GMX_EXPECT_SIMD_REAL_NEAR(refcZ, cZ);
114 }
115
116 #endif      // GMX_SIMD_HAVE_REAL
117
118 /*! \} */
119 /*! \endcond */
120
121 }      // namespace
122 }      // namespace
123 }      // namespace