PME spline+spread CUDA kernel and unit tests
[alexxy/gromacs.git] / src / gromacs / ewald / tests / pmesolvetest.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 /*! \internal \file
36  * \brief
37  * Implements PME solving tests.
38  *
39  * \author Aleksei Iupinov <a.yupinov@gmail.com>
40  * \ingroup module_ewald
41  */
42
43 #include "gmxpre.h"
44
45 #include <string>
46
47 #include <gmock/gmock.h>
48
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/utility/stringutil.h"
51
52 #include "testutils/refdata.h"
53 #include "testutils/testasserts.h"
54
55 #include "pmetestcommon.h"
56
57 namespace gmx
58 {
59 namespace test
60 {
61 namespace
62 {
63 /*! \brief Convenience typedef of the test input parameters - unit cell box, complex grid dimensions, complex grid values,
64  * electrostatic constant epsilon_r, Ewald splitting parameters ewaldcoeff_q and ewaldcoeff_lj, solver type
65  * Output: transformed local grid (Fourier space); optionally reciprocal energy and virial matrix.
66  * TODO:
67  * Implement and test Lorentz-Berthelot
68  */
69 typedef std::tuple<Matrix3x3, IVec, SparseComplexGridValuesInput, double, double, double, PmeSolveAlgorithm> SolveInputParameters;
70
71 //! Test fixture
72 class PmeSolveTest : public ::testing::TestWithParam<SolveInputParameters>
73 {
74     public:
75         //! Default constructor
76         PmeSolveTest() = default;
77
78         //! The test
79         void runTest()
80         {
81             /* Getting the input */
82             Matrix3x3                      box;
83             IVec                           gridSize;
84             SparseComplexGridValuesInput   nonZeroGridValues;
85             double                         epsilon_r;
86             double                         ewaldCoeff_q;
87             double                         ewaldCoeff_lj;
88             PmeSolveAlgorithm              method;
89             std::tie(box, gridSize, nonZeroGridValues, epsilon_r, ewaldCoeff_q, ewaldCoeff_lj, method) = GetParam();
90
91             /* Storing the input where it's needed, running the test */
92             t_inputrec inputRec;
93             inputRec.nkx         = gridSize[XX];
94             inputRec.nky         = gridSize[YY];
95             inputRec.nkz         = gridSize[ZZ];
96             inputRec.pme_order   = 4;
97             inputRec.coulombtype = eelPME;
98             inputRec.epsilon_r   = epsilon_r;
99             switch (method)
100             {
101                 case PmeSolveAlgorithm::Coulomb:
102                     break;
103
104                 case PmeSolveAlgorithm::LennardJones:
105                     inputRec.vdwtype = evdwPME;
106                     break;
107
108                 default:
109                     GMX_THROW(InternalError("Unknown PME solver"));
110             }
111
112             TestReferenceData                     refData;
113             const std::map<CodePath, std::string> modesToTest = {{CodePath::CPU, "CPU"}};
114             for (const auto &mode : modesToTest)
115             {
116                 std::map<GridOrdering, std::string> gridOrderingsToTest = {{GridOrdering::YZX, "YZX"}};
117                 for (const auto &gridOrdering : gridOrderingsToTest)
118                 {
119                     for (bool computeEnergyAndVirial : {false, true})
120                     {
121                         /* Describing the test*/
122                         SCOPED_TRACE(formatString("Testing solving (%s, %s, %s energy/virial) with %s for PME grid size %d %d %d, Ewald coefficients %g %g",
123                                                   (method == PmeSolveAlgorithm::LennardJones) ? "Lennard-Jones" : "Coulomb",
124                                                   gridOrdering.second.c_str(), computeEnergyAndVirial ? "with" : "without",
125                                                   mode.second.c_str(),
126                                                   gridSize[XX], gridSize[YY], gridSize[ZZ],
127                                                   ewaldCoeff_q, ewaldCoeff_lj
128                                                   ));
129
130                         /* Running the test */
131                         PmeSafePointer pmeSafe = pmeInitEmpty(&inputRec, mode.first, nullptr, box, ewaldCoeff_q, ewaldCoeff_lj);
132                         pmeSetComplexGrid(pmeSafe.get(), mode.first, gridOrdering.first, nonZeroGridValues);
133                         const real     cellVolume = box[0] * box[4] * box[8];
134                         //FIXME - this is box[XX][XX] * box[YY][YY] * box[ZZ][ZZ], should be stored in the PME structure
135                         pmePerformSolve(pmeSafe.get(), mode.first, method, cellVolume, gridOrdering.first, computeEnergyAndVirial);
136
137                         /* Check the outputs */
138                         TestReferenceChecker checker(refData.rootChecker());
139                         const auto           ulpTolerance = 200;
140                         checker.setDefaultTolerance(relativeToleranceAsUlp(10.0, ulpTolerance));
141
142                         SparseComplexGridValuesOutput nonZeroGridValuesOutput = pmeGetComplexGrid(pmeSafe.get(), mode.first, gridOrdering.first);
143                         /* Transformed grid */
144                         TestReferenceChecker          gridValuesChecker(checker.checkCompound("NonZeroGridValues", "ComplexSpaceGrid"));
145                         const auto                    ulpToleranceGrid = 50;
146                         gridValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceGrid));
147                         for (const auto &point : nonZeroGridValuesOutput)
148                         {
149                             // we want an additional safeguard for denormal numbers as they cause an exception in string conversion;
150                             // however, using GMX_REAL_MIN causes an "unused item warning" for single precision builds
151                             if (fabs(point.second.re) >= GMX_FLOAT_MIN)
152                             {
153                                 gridValuesChecker.checkReal(point.second.re, (point.first + " re").c_str());
154                             }
155                             if (fabs(point.second.im) >= GMX_FLOAT_MIN)
156                             {
157                                 gridValuesChecker.checkReal(point.second.im, (point.first + " im").c_str());
158                             }
159                         }
160
161                         if (computeEnergyAndVirial)
162                         {
163                             real       energy;
164                             Matrix3x3  virial;
165                             std::tie(energy, virial) = pmeGetReciprocalEnergyAndVirial(pmeSafe.get(), mode.first, method);
166                             /* Energy */
167                             checker.checkReal(energy, "Energy");
168                             /* Virial */
169                             TestReferenceChecker virialChecker(checker.checkCompound("Matrix", "Virial"));
170                             virialChecker.setDefaultTolerance(relativeToleranceAsUlp(1000, 30));
171                             for (int i = 0; i < DIM; i++)
172                             {
173                                 for (int j = 0; j <= i; j++)
174                                 {
175                                     std::string valueId = formatString("Cell %d %d", i, j);
176                                     virialChecker.checkReal(virial[i * DIM + j], valueId.c_str());
177                                 }
178                             }
179                         }
180                     }
181                 }
182             }
183         }
184 };
185
186 /*! \brief Test for PME solving */
187 TEST_P(PmeSolveTest, ReproducesOutputs)
188 {
189     EXPECT_NO_THROW(runTest());
190 }
191
192 /* Valid input instances */
193
194 //! A couple of valid inputs for boxes.
195 static std::vector<Matrix3x3> const c_sampleBoxes
196 {
197     // normal box
198     Matrix3x3 {{
199                    8.0f, 0.0f, 0.0f,
200                    0.0f, 3.4f, 0.0f,
201                    0.0f, 0.0f, 2.0f
202                }},
203     // triclinic box
204     Matrix3x3 {{
205                    7.0f, 0.0f, 0.0f,
206                    0.0f, 4.1f, 0.0f,
207                    3.5f, 2.0f, 12.2f
208                }},
209 };
210
211 //! A couple of valid inputs for grid sizes
212 static std::vector<IVec> const c_sampleGridSizes
213 {
214     IVec {
215         16, 12, 28
216     },
217     IVec {
218         9, 7, 23
219     }
220 };
221
222 //! Moved out from instantiations for readability
223 const auto c_inputBoxes     = ::testing::ValuesIn(c_sampleBoxes);
224 //! Moved out from instantiations for readability
225 const auto c_inputGridSizes = ::testing::ValuesIn(c_sampleGridSizes);
226
227 //! 2 sample complex grids - only non-zero values have to be listed
228 static std::vector<SparseComplexGridValuesInput> const c_sampleGrids
229 {
230     SparseComplexGridValuesInput {{
231                                       IVec {
232                                           0, 0, 0
233                                       }, t_complex {
234                                           3.5f, 6.7f
235                                       }
236                                   }, {
237                                       IVec {
238                                           7, 0, 0
239                                       }, t_complex {
240                                           -2.5f, -0.7f
241                                       }
242                                   }, {
243                                       IVec {
244                                           3, 5, 7
245                                       }, t_complex {
246                                           -0.006f, 1e-8f
247                                       }
248                                   }, {
249                                       IVec {
250                                           3, 1, 2
251                                       }, t_complex {
252                                           0.6f, 7.9f
253                                       }
254                                   },  {
255                                       IVec {
256                                           6, 2, 4
257                                       }, t_complex {
258                                           30.1f, 2.45f
259                                       }
260                                   }, },
261     SparseComplexGridValuesInput {{
262                                       IVec {
263                                           0, 4, 0
264                                       }, t_complex {
265                                           0.0f, 0.3f
266                                       }
267                                   }, {
268                                       IVec {
269                                           4, 2, 7
270                                       }, t_complex {
271                                           13.76f, -40.0f
272                                       }
273                                   }, {
274                                       IVec {
275                                           0, 6, 7
276                                       }, t_complex {
277                                           3.6f, 0.0f
278                                       }
279                                   }, {
280                                       IVec {
281                                           2, 5, 10
282                                       }, t_complex {
283                                           3.6f, 10.65f
284                                       }
285                                   }, }
286 };
287
288 //! Moved out from instantiations for readability
289 const auto c_inputGrids = ::testing::ValuesIn(c_sampleGrids);
290 //! Moved out from instantiations for readability
291 const auto c_inputEpsilon_r = ::testing::Values(1.2);
292 //! Moved out from instantiations for readability
293 const auto c_inputEwaldCoeff_q = ::testing::Values(2.0);
294 //! Moved out from instantiations for readability
295 const auto c_inputEwaldCoeff_lj = ::testing::Values(0.7);
296 //! Moved out from instantiations for readability
297 const auto c_inputMethods = ::testing::Values(PmeSolveAlgorithm::Coulomb, PmeSolveAlgorithm::LennardJones);
298
299 //! Instantiation of the PME solving test
300 INSTANTIATE_TEST_CASE_P(SaneInput, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
301                                                                     c_inputEpsilon_r, c_inputEwaldCoeff_q, c_inputEwaldCoeff_lj, c_inputMethods));
302
303 //! A few more instances to check that different ewaldCoeff_q actually affects results of the Coulomb solver
304 INSTANTIATE_TEST_CASE_P(DifferentEwaldCoeffQ, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
305                                                                                c_inputEpsilon_r, ::testing::Values(0.4), c_inputEwaldCoeff_lj,
306                                                                                    ::testing::Values(PmeSolveAlgorithm::Coulomb)));
307
308 //! A few more instances to check that different ewaldCoeff_lj actually affects results of the Lennard-Jones solver.
309 //! The value has to be approximately larger than 1 / (box dimensions) to have a meaningful output grid.
310 //! Previous value of 0.3 caused one of the grid cells to be less or greater than GMX_FLOAT_MIN, depending on the architecture.
311 INSTANTIATE_TEST_CASE_P(DifferentEwaldCoeffLJ, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
312                                                                                 c_inputEpsilon_r, c_inputEwaldCoeff_q, ::testing::Values(2.35),
313                                                                                     ::testing::Values(PmeSolveAlgorithm::LennardJones)));
314
315 //! A few more instances to check that different epsilon_r actually affects results of all solvers
316 INSTANTIATE_TEST_CASE_P(DifferentEpsilonR, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
317                                                                             testing::Values(1.9), c_inputEwaldCoeff_q, c_inputEwaldCoeff_lj,
318                                                                             c_inputMethods));
319
320 }
321 }
322 }