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