Merge branch release-2016 into release-2018
[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                                                                  {CodePath::CUDA, "CUDA"}};
115             for (const auto &mode : modesToTest)
116             {
117                 const bool supportedInput = pmeSupportsInputForMode(&inputRec, mode.first);
118                 if (!supportedInput)
119                 {
120                     /* Testing the failure for the unsupported input */
121                     EXPECT_THROW(pmeInitEmpty(&inputRec, mode.first, nullptr, box, ewaldCoeff_q, ewaldCoeff_lj), NotImplementedError);
122                     continue;
123                 }
124
125                 std::map<GridOrdering, std::string> gridOrderingsToTest = {{GridOrdering::YZX, "YZX"}};
126                 if (mode.first == CodePath::CUDA)
127                 {
128                     gridOrderingsToTest[GridOrdering::XYZ] = "XYZ";
129                 }
130                 const auto contextsToTest = getPmeTestEnv()->getHardwareContexts(mode.first);
131                 for (const auto &gridOrdering : gridOrderingsToTest)
132                 {
133                     for (const auto &context : contextsToTest)
134                     {
135                         for (bool computeEnergyAndVirial : {false, true})
136                         {
137                             /* Describing the test*/
138                             SCOPED_TRACE(formatString("Testing solving (%s, %s, %s energy/virial) with %s %sfor PME grid size %d %d %d, Ewald coefficients %g %g",
139                                                       (method == PmeSolveAlgorithm::LennardJones) ? "Lennard-Jones" : "Coulomb",
140                                                       gridOrdering.second.c_str(),
141                                                       computeEnergyAndVirial ? "with" : "without",
142                                                       mode.second.c_str(),
143                                                       context.getDescription().c_str(),
144                                                       gridSize[XX], gridSize[YY], gridSize[ZZ],
145                                                       ewaldCoeff_q, ewaldCoeff_lj
146                                                       ));
147
148                             /* Running the test */
149                             PmeSafePointer pmeSafe = pmeInitEmpty(&inputRec, mode.first, context.getDeviceInfo(), box, ewaldCoeff_q, ewaldCoeff_lj);
150                             pmeSetComplexGrid(pmeSafe.get(), mode.first, gridOrdering.first, nonZeroGridValues);
151                             const real     cellVolume = box[0] * box[4] * box[8];
152                             //FIXME - this is box[XX][XX] * box[YY][YY] * box[ZZ][ZZ], should be stored in the PME structure
153                             pmePerformSolve(pmeSafe.get(), mode.first, method, cellVolume, gridOrdering.first, computeEnergyAndVirial);
154                             pmeFinalizeTest(pmeSafe.get(), mode.first);
155
156                             /* Check the outputs */
157                             TestReferenceChecker          checker(refData.rootChecker());
158
159                             SparseComplexGridValuesOutput nonZeroGridValuesOutput = pmeGetComplexGrid(pmeSafe.get(), mode.first, gridOrdering.first);
160                             /* Transformed grid */
161                             TestReferenceChecker          gridValuesChecker(checker.checkCompound("NonZeroGridValues", "ComplexSpaceGrid"));
162                             const auto                    ulpToleranceGrid = 40;
163                             gridValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceGrid));
164                             for (const auto &point : nonZeroGridValuesOutput)
165                             {
166                                 // we want an additional safeguard for denormal numbers as they cause an exception in string conversion;
167                                 // however, using GMX_REAL_MIN causes an "unused item warning" for single precision builds
168                                 if (fabs(point.second.re) >= GMX_FLOAT_MIN)
169                                 {
170                                     gridValuesChecker.checkReal(point.second.re, (point.first + " re").c_str());
171                                 }
172                                 if (fabs(point.second.im) >= GMX_FLOAT_MIN)
173                                 {
174                                     gridValuesChecker.checkReal(point.second.im, (point.first + " im").c_str());
175                                 }
176                             }
177
178                             if (computeEnergyAndVirial)
179                             {
180                                 real       energy;
181                                 Matrix3x3  virial;
182                                 std::tie(energy, virial) = pmeGetReciprocalEnergyAndVirial(pmeSafe.get(), mode.first, method);
183                                 /* Energy */
184                                 TestReferenceChecker energyChecker(checker);
185                                 energyChecker.setDefaultTolerance(relativeToleranceAsUlp(10.0, 40));
186                                 energyChecker.checkReal(energy, "Energy");
187                                 /* Virial */
188                                 TestReferenceChecker virialChecker(checker.checkCompound("Matrix", "Virial"));
189                                 virialChecker.setDefaultTolerance(relativeToleranceAsUlp(1000, 8));
190                                 for (int i = 0; i < DIM; i++)
191                                 {
192                                     for (int j = 0; j <= i; j++)
193                                     {
194                                         std::string valueId = formatString("Cell %d %d", i, j);
195                                         virialChecker.checkReal(virial[i * DIM + j], valueId.c_str());
196                                     }
197                                 }
198                             }
199                         }
200                     }
201                 }
202
203             }
204         }
205 };
206
207 /*! \brief Test for PME solving */
208 TEST_P(PmeSolveTest, ReproducesOutputs)
209 {
210     EXPECT_NO_THROW(runTest());
211 }
212
213 /* Valid input instances */
214
215 //! A couple of valid inputs for boxes.
216 static std::vector<Matrix3x3> const c_sampleBoxes
217 {
218     // normal box
219     Matrix3x3 {{
220                    8.0f, 0.0f, 0.0f,
221                    0.0f, 3.4f, 0.0f,
222                    0.0f, 0.0f, 2.0f
223                }},
224     // triclinic box
225     Matrix3x3 {{
226                    7.0f, 0.0f, 0.0f,
227                    0.0f, 4.1f, 0.0f,
228                    3.5f, 2.0f, 12.2f
229                }},
230 };
231
232 //! A couple of valid inputs for grid sizes
233 static std::vector<IVec> const c_sampleGridSizes
234 {
235     IVec {
236         16, 12, 28
237     },
238     IVec {
239         9, 7, 23
240     }
241 };
242
243 //! Moved out from instantiations for readability
244 const auto c_inputBoxes     = ::testing::ValuesIn(c_sampleBoxes);
245 //! Moved out from instantiations for readability
246 const auto c_inputGridSizes = ::testing::ValuesIn(c_sampleGridSizes);
247
248 //! 2 sample complex grids - only non-zero values have to be listed
249 static std::vector<SparseComplexGridValuesInput> const c_sampleGrids
250 {
251     SparseComplexGridValuesInput {{
252                                       IVec {
253                                           0, 0, 0
254                                       }, t_complex {
255                                           3.5f, 6.7f
256                                       }
257                                   }, {
258                                       IVec {
259                                           7, 0, 0
260                                       }, t_complex {
261                                           -2.5f, -0.7f
262                                       }
263                                   }, {
264                                       IVec {
265                                           3, 5, 7
266                                       }, t_complex {
267                                           -0.006f, 1e-8f
268                                       }
269                                   }, {
270                                       IVec {
271                                           3, 1, 2
272                                       }, t_complex {
273                                           0.6f, 7.9f
274                                       }
275                                   },  {
276                                       IVec {
277                                           6, 2, 4
278                                       }, t_complex {
279                                           30.1f, 2.45f
280                                       }
281                                   }, },
282     SparseComplexGridValuesInput {{
283                                       IVec {
284                                           0, 4, 0
285                                       }, t_complex {
286                                           0.0f, 0.3f
287                                       }
288                                   }, {
289                                       IVec {
290                                           4, 2, 7
291                                       }, t_complex {
292                                           13.76f, -40.0f
293                                       }
294                                   }, {
295                                       IVec {
296                                           0, 6, 7
297                                       }, t_complex {
298                                           3.6f, 0.0f
299                                       }
300                                   }, {
301                                       IVec {
302                                           2, 5, 10
303                                       }, t_complex {
304                                           3.6f, 10.65f
305                                       }
306                                   }, }
307 };
308
309 //! Moved out from instantiations for readability
310 const auto c_inputGrids = ::testing::ValuesIn(c_sampleGrids);
311 //! Moved out from instantiations for readability
312 const auto c_inputEpsilon_r = ::testing::Values(1.2);
313 //! Moved out from instantiations for readability
314 const auto c_inputEwaldCoeff_q = ::testing::Values(2.0);
315 //! Moved out from instantiations for readability
316 const auto c_inputEwaldCoeff_lj = ::testing::Values(0.7);
317 //! Moved out from instantiations for readability
318 const auto c_inputMethods = ::testing::Values(PmeSolveAlgorithm::Coulomb, PmeSolveAlgorithm::LennardJones);
319
320 //! Instantiation of the PME solving test
321 INSTANTIATE_TEST_CASE_P(SaneInput, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
322                                                                     c_inputEpsilon_r, c_inputEwaldCoeff_q, c_inputEwaldCoeff_lj, c_inputMethods));
323
324 //! A few more instances to check that different ewaldCoeff_q actually affects results of the Coulomb solver
325 INSTANTIATE_TEST_CASE_P(DifferentEwaldCoeffQ, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
326                                                                                c_inputEpsilon_r, ::testing::Values(0.4), c_inputEwaldCoeff_lj,
327                                                                                    ::testing::Values(PmeSolveAlgorithm::Coulomb)));
328
329 //! A few more instances to check that different ewaldCoeff_lj actually affects results of the Lennard-Jones solver.
330 //! The value has to be approximately larger than 1 / (box dimensions) to have a meaningful output grid.
331 //! Previous value of 0.3 caused one of the grid cells to be less or greater than GMX_FLOAT_MIN, depending on the architecture.
332 INSTANTIATE_TEST_CASE_P(DifferentEwaldCoeffLJ, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
333                                                                                 c_inputEpsilon_r, c_inputEwaldCoeff_q, ::testing::Values(2.35),
334                                                                                     ::testing::Values(PmeSolveAlgorithm::LennardJones)));
335
336 //! A few more instances to check that different epsilon_r actually affects results of all solvers
337 INSTANTIATE_TEST_CASE_P(DifferentEpsilonR, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
338                                                                             testing::Values(1.9), c_inputEwaldCoeff_q, c_inputEwaldCoeff_lj,
339                                                                             c_inputMethods));
340
341 }
342 }
343 }