Extend support for float/double tolerances in testing
[alexxy/gromacs.git] / src / gromacs / fft / tests / fft.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,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 /*! \internal \file
36  * \brief
37  * Tests utilities for fft calculations.
38  *
39  * Current reference data is generated in double precision using the Reference
40  * build type, except for the compiler (Apple Clang).
41  *
42  * \author Roland Schulz <roland@utk.edu>
43  * \ingroup module_fft
44  */
45 #include "gmxpre.h"
46
47 #include "gromacs/fft/fft.h"
48
49 #include <vector>
50
51 #include <gtest/gtest.h>
52
53 #include "gromacs/fft/parallel_3dfft.h"
54 #include "gromacs/utility/stringutil.h"
55
56 #include "testutils/refdata.h"
57 #include "testutils/testasserts.h"
58
59 namespace
60 {
61
62 //! Input data for FFT tests.
63 const real inputdata[] = { //print ",\n".join([",".join(["%4s"%(random.randint(-99,99)/10.,) for i in range(25)]) for j in range(20)])
64     -3.5, 6.3, 1.2, 0.3, 1.1, -5.7, 5.8, -1.9, -6.3, -1.4, 7.4, 2.4, -9.9, -7.2, 5.4, 6.1, -1.9, -7.6, 1.4, -3.5, 0.7, 5.6, -4.2, -1.1, -4.4,
65     -6.3, -7.2, 4.6, -3.0, -0.9, 7.2, 2.5, -3.6, 6.1, -3.2, -2.1, 6.5, -0.4, -9.0, 2.3, 8.4, 4.0, -5.2, -9.0, 4.7, -3.7, -2.0, -9.5, -3.9, -3.6,
66     7.1, 0.8, -0.6, 5.2, -9.3, -4.5, 5.9, 2.2, -5.8, 5.0, 1.2, -0.1, 2.2, 0.2, -7.7, 1.9, -8.4, 4.4, 2.3, -2.9, 6.7, 2.7, 5.8, -3.6, 8.9,
67     8.9, 4.3, 9.1, 9.3, -8.7, 4.1, 9.6, -6.2, 6.6, -9.3, 8.2, 4.5, 6.2, 9.4, -8.0, -6.8, -3.3, 7.2, 1.7, 0.6, -4.9, 9.8, 1.3, 3.2, -0.2,
68     9.9, 4.4, -9.9, -7.2, 4.4, 4.7, 7.2, -0.3, 0.3, -2.1, 8.4, -2.1, -6.1, 4.1, -5.9, -2.2, -3.8, 5.2, -8.2, -7.8, -8.8, 6.7, -9.5, -4.2, 0.8,
69     8.3, 5.2, -9.0, 8.7, 9.8, -9.9, -7.8, -8.3, 9.0, -2.8, -9.2, -9.6, 8.4, 2.5, 6.0, -0.4, 1.3, -0.5, 9.1, -9.5, -0.8, 1.9, -6.2, 4.3, -3.8,
70     8.6, -1.9, -2.1, -0.4, -7.1, -3.7, 9.1, -6.4, -0.6, 2.5, 8.0, -5.2, -9.8, -4.3, 4.5, 1.7, 9.3, 9.2, 1.0, 5.3, -4.5, 6.4, -6.6, 3.1, -6.8,
71     2.1, 2.0, 7.3, 8.6, 5.0, 5.2, 0.4, -7.1, 4.5, -9.2, -9.1, 0.2, -6.3, -1.1, -9.6, 7.4, -3.7, -5.5, 2.6, -3.5, -0.7, 9.0, 9.8, -8.0, 3.6,
72     3.0, -2.2, -2.8, 0.8, 9.0, 2.8, 7.7, -0.7, -5.0, -1.8, -2.3, -0.4, -6.2, -9.1, -9.2, 0.5, 5.7, -3.9, 2.1, 0.6, 0.4, 9.1, 7.4, 7.1, -2.5,
73     7.3, 7.8, -4.3, 6.3, -0.8, -3.8, -1.5, 6.6, 2.3, 3.9, -4.6, 5.8, -7.4, 5.9, 2.8, 4.7, 3.9, -5.4, 9.1, -1.6, -1.9, -4.2, -2.6, 0.6, -5.1,
74     1.8, 5.2, 4.0, -6.2, 6.5, -9.1, 0.5, 2.1, 7.1, -8.6, 7.6, -9.7, -4.6, -5.7, 6.1, -1.8, -7.3, 9.4, 8.0, -2.6, -1.8, 5.7, 9.3, -7.9, 7.4,
75     6.3, 2.0, 9.6, -4.5, -6.2, 6.1, 2.3, 0.8, 5.9, -2.8, -3.5, -1.5, 6.0, -4.9, 3.5, 7.7, -4.2, -9.7, 2.4, 8.1, 5.9, 3.4, -7.5, 7.5, 2.6,
76     4.7, 2.7, 2.2, 2.6, 6.2, 7.5, 0.2, -6.4, -2.8, -0.5, -0.3, 0.4, 1.2, 3.5, -4.0, -0.5, 9.3, -7.2, 8.5, -5.5, -1.7, -5.3, 0.3, 3.9, -3.6,
77     -3.6, 4.7, -8.1, 1.4, 4.0, 1.3, -4.3, -8.8, -7.3, 6.3, -7.5, -9.0, 9.1, 4.5, -1.9, 1.9, 9.9, -1.7, -9.1, -5.1, 8.5, -9.3, 2.1, -5.8, -3.6,
78     -0.8, -0.9, -3.3, -2.7, 7.0, -7.2, -5.0, 7.4, -1.4, 0.0, -4.5, -9.7, 0.7, -1.0, -9.1, -5.3, 4.3, 3.4, -6.6, 9.8, -1.1, 8.9, 5.0, 2.9, 0.2,
79     -2.9, 0.8, 6.7, -0.6, 0.6, 4.1, 5.3, -1.7, -0.3, 4.2, 3.7, -8.3, 4.0, 1.3, 6.3, 0.2, 1.3, -1.1, -3.5, 2.8, -7.7, 6.2, -4.9, -9.9, 9.6,
80     3.0, -9.2, -8.0, -3.9, 7.9, -6.1, 6.0, 5.9, 9.6, 1.2, 6.2, 3.6, 2.1, 5.8, 9.2, -8.8, 8.8, -3.3, -9.2, 4.6, 1.8, 4.6, 2.9, -2.7, 4.2,
81     7.3, -0.4, 7.7, -7.0, 2.1, 0.3, 3.7, 3.3, -8.6, 9.8, 3.6, 3.1, 6.5, -2.4, 7.8, 7.5, 8.4, -2.8, -6.3, -5.1, -2.7, 9.3, -0.8, -9.2, 7.9,
82     8.9, 3.4, 0.1, -5.3, -6.8, 4.9, 4.3, -0.7, -2.2, -3.2, -7.5, -2.3, 0.0, 8.1, -9.2, -2.3, -5.7, 2.1, 2.6, 2.0, 0.3, -8.0, -2.0, -7.9, 6.6,
83     8.4, 4.0, -6.2, -6.9, -7.2, 7.7, -5.0, 5.3, 1.9, -5.3, -7.5, 8.8, 8.3, 9.0, 8.1, 3.2, 1.2, -5.4, -0.2, 2.1, -5.2, 9.5, 5.9, 5.6, -7.8,
84 };
85
86
87 class BaseFFTTest : public ::testing::Test
88 {
89     public:
90         BaseFFTTest()
91             : checker_(data_.rootChecker()), flags_(GMX_FFT_FLAG_CONSERVATIVE)
92         {
93             // TODO: These tolerances are just something that has been observed
94             // to be sufficient to pass the tests.  It would be nicer to
95             // actually argue about why they are sufficient (or what is).
96             checker_.setDefaultTolerance(
97                     gmx::test::relativeToleranceAsPrecisionDependentUlp(10.0, 64, 512));
98         }
99         ~BaseFFTTest()
100         {
101             gmx_fft_cleanup();
102         }
103
104         gmx::test::TestReferenceData    data_;
105         gmx::test::TestReferenceChecker checker_;
106         std::vector<real>               in_, out_;
107         int                             flags_;
108 };
109
110 class FFTTest : public BaseFFTTest
111 {
112     public:
113         FFTTest() : fft_(NULL)
114         {
115         }
116         ~FFTTest()
117         {
118             if (fft_)
119             {
120                 gmx_fft_destroy(fft_);
121             }
122         }
123         gmx_fft_t fft_;
124 };
125
126 class ManyFFTTest : public BaseFFTTest
127 {
128     public:
129         ManyFFTTest() : fft_(NULL)
130         {
131         }
132         ~ManyFFTTest()
133         {
134             if (fft_)
135             {
136                 gmx_many_fft_destroy(fft_);
137             }
138         }
139         gmx_fft_t fft_;
140 };
141
142
143 //TODO: Add tests for aligned/not-aligned input/output memory
144
145 class FFTTest1D : public FFTTest, public ::testing::WithParamInterface<int>
146 {
147
148 };
149
150 class FFFTest3D : public BaseFFTTest
151 {
152     public:
153         FFFTest3D() : fft_(NULL)
154         {
155         }
156         ~FFFTest3D()
157         {
158             if (fft_)
159             {
160                 gmx_parallel_3dfft_destroy(fft_);
161             }
162         }
163         gmx_parallel_3dfft_t fft_;
164 };
165
166
167 TEST_P(FFTTest1D, Complex)
168 {
169     const int nx = GetParam();
170     ASSERT_LE(nx*2, static_cast<int>(sizeof(inputdata)/sizeof(real)));
171
172     in_  = std::vector<real>(inputdata, inputdata+nx*2);
173     out_ = std::vector<real>(nx*2);
174     real* in  = &in_[0];
175     real* out = &out_[0];
176
177     gmx_fft_init_1d(&fft_, nx, flags_);
178
179     gmx_fft_1d(fft_, GMX_FFT_FORWARD, in, out);
180     checker_.checkSequenceArray(nx*2, out, "forward");
181     gmx_fft_1d(fft_, GMX_FFT_BACKWARD, in, out);
182     checker_.checkSequenceArray(nx*2, out, "backward");
183 }
184
185 TEST_P(FFTTest1D, Real)
186 {
187     const int rx = GetParam();
188     const int cx = (rx/2+1);
189     ASSERT_LE(cx*2, static_cast<int>(sizeof(inputdata)/sizeof(real)));
190
191     in_  = std::vector<real>(inputdata, inputdata+cx*2);
192     out_ = std::vector<real>(cx*2);
193     real* in  = &in_[0];
194     real* out = &out_[0];
195
196     gmx_fft_init_1d_real(&fft_, rx, flags_);
197
198     gmx_fft_1d_real(fft_, GMX_FFT_REAL_TO_COMPLEX, in, out);
199     checker_.checkSequenceArray(cx*2, out, "forward");
200     gmx_fft_1d_real(fft_, GMX_FFT_COMPLEX_TO_REAL, in, out);
201     checker_.checkSequenceArray(rx, out, "backward");
202 }
203
204 INSTANTIATE_TEST_CASE_P(7_8_25_36_60,
205                         FFTTest1D, ::testing::Values(7, 8, 25, 36, 60));
206
207
208 TEST_F(ManyFFTTest, Complex1DLength48Multi5Test)
209 {
210     const int nx = 48;
211     const int N  = 5;
212
213     in_  = std::vector<real>(inputdata, inputdata+nx*2*N);
214     out_ = std::vector<real>(nx*2*N);
215     real* in  = &in_[0];
216     real* out = &out_[0];
217
218     gmx_fft_init_many_1d(&fft_, nx, N, flags_);
219
220     gmx_fft_many_1d(fft_, GMX_FFT_FORWARD, in, out);
221     checker_.checkSequenceArray(nx*2*N, out, "forward");
222     gmx_fft_many_1d(fft_, GMX_FFT_BACKWARD, in, out);
223     checker_.checkSequenceArray(nx*2*N, out, "backward");
224 }
225
226 TEST_F(ManyFFTTest, Real1DLength48Multi5Test)
227 {
228     const int rx = 48;
229     const int cx = (rx/2+1);
230     const int N  = 5;
231
232     in_  = std::vector<real>(inputdata, inputdata+cx*2*N);
233     out_ = std::vector<real>(cx*2*N);
234     real* in  = &in_[0];
235     real* out = &out_[0];
236
237     gmx_fft_init_many_1d_real(&fft_, rx, N, flags_);
238
239     gmx_fft_many_1d_real(fft_, GMX_FFT_REAL_TO_COMPLEX, in, out);
240     checker_.checkSequenceArray(cx*2*N, out, "forward");
241     gmx_fft_many_1d_real(fft_, GMX_FFT_COMPLEX_TO_REAL, in, out);
242     checker_.checkSequenceArray(rx*N, out, "backward");
243 }
244
245 TEST_F(FFTTest, Real2DLength18_15Test)
246 {
247     const int rx = 18;
248     const int cx = (rx/2+1);
249     const int ny = 15;
250
251     in_  = std::vector<real>(inputdata, inputdata+cx*2*ny);
252     out_ = std::vector<real>(cx*2*ny);
253     real* in  = &in_[0];
254     real* out = &out_[0];
255
256     gmx_fft_init_2d_real(&fft_, rx, ny, flags_);
257
258     gmx_fft_2d_real(fft_, GMX_FFT_REAL_TO_COMPLEX, in, out);
259     checker_.checkSequenceArray(cx*2*ny, out, "forward");
260 //    known to be wrong for gmx_fft_mkl. And not used.
261 //    gmx_fft_2d_real(_fft,GMX_FFT_COMPLEX_TO_REAL,in,out);
262 //    _checker.checkSequenceArray(rx*ny, out, "backward");
263 }
264
265 //TODO: test with threads and more than 1 MPI ranks
266 TEST_F(FFFTest3D, Real5_6_9)
267 {
268     int        ndata[] = {5, 6, 9};
269     MPI_Comm   comm[]  = {MPI_COMM_NULL, MPI_COMM_NULL};
270     real     * rdata;
271     t_complex* cdata;
272     ivec       local_ndata, offset, rsize, csize, complex_order;
273
274     gmx_parallel_3dfft_init(&fft_, ndata, &rdata, &cdata,
275                             comm, TRUE, 1);
276
277     gmx_parallel_3dfft_real_limits(fft_, local_ndata, offset, rsize);
278     gmx_parallel_3dfft_complex_limits(fft_, complex_order,
279                                       local_ndata, offset, csize);
280     checker_.checkVector(rsize, "rsize");
281     checker_.checkVector(csize, "csize");
282     int size = csize[0]*csize[1]*csize[2];
283
284     memcpy(rdata, inputdata, size*sizeof(t_complex));
285     gmx_parallel_3dfft_execute(fft_, GMX_FFT_REAL_TO_COMPLEX, 0, NULL);
286     //TODO use std::complex and add checkComplex for it
287     checker_.checkSequenceArray(size*2,
288                                 reinterpret_cast<real*>(cdata), "forward");
289
290     memcpy(cdata, inputdata, size*sizeof(t_complex));
291     gmx_parallel_3dfft_execute(fft_, GMX_FFT_COMPLEX_TO_REAL, 0, NULL);
292     for (int i = 0; i < ndata[0]*ndata[1]; i++) //check sequence but skip unused data
293     {
294         checker_.checkSequenceArray(ndata[2], rdata+i*rsize[2],
295                                     gmx::formatString("backward %d", i).c_str());
296     }
297 }
298
299 } // namespace