5c48955cc473f27f37737dce24a47f83052b0eb3
[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,2019,2020,2021, 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/test_hardware_environment.h"
54 #include "testutils/testasserts.h"
55
56 #include "pmetestcommon.h"
57
58 namespace gmx
59 {
60 namespace test
61 {
62 namespace
63 {
64 /*! \brief Convenience typedef of the test input parameters - unit cell box, complex grid dimensions, complex grid values,
65  * electrostatic constant epsilon_r, Ewald splitting parameters ewaldcoeff_q and ewaldcoeff_lj, solver type
66  * Output: transformed local grid (Fourier space); optionally reciprocal energy and virial matrix.
67  * TODO:
68  * Implement and test Lorentz-Berthelot
69  */
70 typedef std::tuple<Matrix3x3, IVec, SparseComplexGridValuesInput, double, double, double, PmeSolveAlgorithm> SolveInputParameters;
71
72 //! Test fixture
73 class PmeSolveTest : public ::testing::TestWithParam<SolveInputParameters>
74 {
75 public:
76     PmeSolveTest() = default;
77
78     //! Sets the programs once
79     static void SetUpTestSuite() { s_pmeTestHardwareContexts = createPmeTestHardwareContextList(); }
80
81     //! The test
82     static void runTest()
83     {
84         /* Getting the input */
85         Matrix3x3                    box;
86         IVec                         gridSize;
87         SparseComplexGridValuesInput nonZeroGridValues;
88         double                       epsilon_r;
89         double                       ewaldCoeff_q;
90         double                       ewaldCoeff_lj;
91         PmeSolveAlgorithm            method;
92         std::tie(box, gridSize, nonZeroGridValues, epsilon_r, ewaldCoeff_q, ewaldCoeff_lj, method) =
93                 GetParam();
94
95         /* Storing the input where it's needed, running the test */
96         t_inputrec inputRec;
97         inputRec.nkx         = gridSize[XX];
98         inputRec.nky         = gridSize[YY];
99         inputRec.nkz         = gridSize[ZZ];
100         inputRec.pme_order   = 4;
101         inputRec.coulombtype = CoulombInteractionType::Pme;
102         inputRec.epsilon_r   = epsilon_r;
103         switch (method)
104         {
105             case PmeSolveAlgorithm::Coulomb: break;
106
107             case PmeSolveAlgorithm::LennardJones: inputRec.vdwtype = VanDerWaalsType::Pme; break;
108
109             default: GMX_THROW(InternalError("Unknown PME solver"));
110         }
111
112         TestReferenceData refData;
113         for (const auto& pmeTestHardwareContext : s_pmeTestHardwareContexts)
114         {
115             pmeTestHardwareContext->activate();
116             CodePath   codePath       = pmeTestHardwareContext->codePath();
117             const bool supportedInput = pmeSupportsInputForMode(
118                     *getTestHardwareEnvironment()->hwinfo(), &inputRec, codePath);
119             if (!supportedInput)
120             {
121                 /* Testing the failure for the unsupported input */
122                 EXPECT_THROW_GMX(
123                         pmeInitWrapper(&inputRec, codePath, nullptr, nullptr, nullptr, box, ewaldCoeff_q, ewaldCoeff_lj),
124                         NotImplementedError);
125                 continue;
126             }
127
128             std::map<GridOrdering, std::string> gridOrderingsToTest = { { GridOrdering::YZX,
129                                                                           "YZX" } };
130             if (codePath == CodePath::GPU)
131             {
132                 gridOrderingsToTest[GridOrdering::XYZ] = "XYZ";
133             }
134             for (const auto& gridOrdering : gridOrderingsToTest)
135             {
136                 for (bool computeEnergyAndVirial : { false, true })
137                 {
138                     /* Describing the test*/
139                     SCOPED_TRACE(formatString(
140                             "Testing solving (%s, %s, %s energy/virial) on %s for PME grid "
141                             "size %d %d %d, Ewald coefficients %g %g",
142                             (method == PmeSolveAlgorithm::LennardJones) ? "Lennard-Jones" : "Coulomb",
143                             gridOrdering.second.c_str(),
144                             computeEnergyAndVirial ? "with" : "without",
145                             pmeTestHardwareContext->description().c_str(),
146                             gridSize[XX],
147                             gridSize[YY],
148                             gridSize[ZZ],
149                             ewaldCoeff_q,
150                             ewaldCoeff_lj));
151
152                     /* Running the test */
153                     PmeSafePointer pmeSafe = pmeInitWrapper(&inputRec,
154                                                             codePath,
155                                                             pmeTestHardwareContext->deviceContext(),
156                                                             pmeTestHardwareContext->deviceStream(),
157                                                             pmeTestHardwareContext->pmeGpuProgram(),
158                                                             box,
159                                                             ewaldCoeff_q,
160                                                             ewaldCoeff_lj);
161                     pmeSetComplexGrid(pmeSafe.get(), codePath, gridOrdering.first, nonZeroGridValues);
162                     const real cellVolume = box[0] * box[4] * box[8];
163                     // FIXME - this is box[XX][XX] * box[YY][YY] * box[ZZ][ZZ], should be stored in the PME structure
164                     pmePerformSolve(pmeSafe.get(), codePath, method, cellVolume, gridOrdering.first, computeEnergyAndVirial);
165                     pmeFinalizeTest(pmeSafe.get(), codePath);
166
167                     /* Check the outputs */
168                     TestReferenceChecker checker(refData.rootChecker());
169
170                     SparseComplexGridValuesOutput nonZeroGridValuesOutput =
171                             pmeGetComplexGrid(pmeSafe.get(), codePath, gridOrdering.first);
172                     /* Transformed grid */
173                     TestReferenceChecker gridValuesChecker(
174                             checker.checkCompound("NonZeroGridValues", "ComplexSpaceGrid"));
175
176                     real gridValuesMagnitude = 1.0;
177                     for (const auto& point : nonZeroGridValuesOutput)
178                     {
179                         gridValuesMagnitude = std::max(std::fabs(point.second.re), gridValuesMagnitude);
180                         gridValuesMagnitude = std::max(std::fabs(point.second.im), gridValuesMagnitude);
181                     }
182                     // Spline moduli participate 3 times in the computation; 2 is an additional factor for SIMD exp() precision
183                     uint64_t gridUlpToleranceFactor = DIM * 2;
184                     if (method == PmeSolveAlgorithm::LennardJones)
185                     {
186                         // Lennard Jones is more complex and also uses erfc(), relax more
187                         gridUlpToleranceFactor *= 2;
188                     }
189                     const uint64_t splineModuliDoublePrecisionUlps =
190                             getSplineModuliDoublePrecisionUlps(inputRec.pme_order + 1);
191                     auto gridTolerance = relativeToleranceAsPrecisionDependentUlp(
192                             gridValuesMagnitude,
193                             gridUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
194                             gridUlpToleranceFactor * splineModuliDoublePrecisionUlps);
195                     gridValuesChecker.setDefaultTolerance(gridTolerance);
196
197                     for (const auto& point : nonZeroGridValuesOutput)
198                     {
199                         // we want an additional safeguard for denormal numbers as they cause an exception in string conversion;
200                         // however, using GMX_REAL_MIN causes an "unused item warning" for single precision builds
201                         if (fabs(point.second.re) >= GMX_FLOAT_MIN)
202                         {
203                             gridValuesChecker.checkReal(point.second.re, (point.first + " re").c_str());
204                         }
205                         if (fabs(point.second.im) >= GMX_FLOAT_MIN)
206                         {
207                             gridValuesChecker.checkReal(point.second.im, (point.first + " im").c_str());
208                         }
209                     }
210
211                     if (computeEnergyAndVirial)
212                     {
213                         // Extract the energy and virial
214                         const auto output =
215                                 pmeGetReciprocalEnergyAndVirial(pmeSafe.get(), codePath, method);
216                         const auto& energy = (method == PmeSolveAlgorithm::Coulomb)
217                                                      ? output.coulombEnergy_
218                                                      : output.lennardJonesEnergy_;
219                         const auto& virial = (method == PmeSolveAlgorithm::Coulomb)
220                                                      ? output.coulombVirial_
221                                                      : output.lennardJonesVirial_;
222
223                         // These quantities are computed based on the grid values, so must have
224                         // checking relative tolerances at least as large. Virial needs more flops
225                         // than energy, so needs a larger tolerance.
226
227                         /* Energy */
228                         double energyMagnitude = 10.0;
229                         // TODO This factor is arbitrary, do a proper error-propagation analysis
230                         uint64_t energyUlpToleranceFactor = gridUlpToleranceFactor * 2;
231                         auto     energyTolerance = relativeToleranceAsPrecisionDependentUlp(
232                                 energyMagnitude,
233                                 energyUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
234                                 energyUlpToleranceFactor * splineModuliDoublePrecisionUlps);
235                         TestReferenceChecker energyChecker(checker);
236                         energyChecker.setDefaultTolerance(energyTolerance);
237                         energyChecker.checkReal(energy, "Energy");
238
239                         /* Virial */
240                         double virialMagnitude = 1000.0;
241                         // TODO This factor is arbitrary, do a proper error-propagation analysis
242                         uint64_t virialUlpToleranceFactor = energyUlpToleranceFactor * 2;
243                         auto     virialTolerance = relativeToleranceAsPrecisionDependentUlp(
244                                 virialMagnitude,
245                                 virialUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
246                                 virialUlpToleranceFactor * splineModuliDoublePrecisionUlps);
247                         TestReferenceChecker virialChecker(
248                                 checker.checkCompound("Matrix", "Virial"));
249                         virialChecker.setDefaultTolerance(virialTolerance);
250                         for (int i = 0; i < DIM; i++)
251                         {
252                             for (int j = 0; j <= i; j++)
253                             {
254                                 std::string valueId = formatString("Cell %d %d", i, j);
255                                 virialChecker.checkReal(virial[i][j], valueId.c_str());
256                             }
257                         }
258                     }
259                 }
260             }
261         }
262     }
263
264     static std::vector<std::unique_ptr<PmeTestHardwareContext>> s_pmeTestHardwareContexts;
265 };
266
267 std::vector<std::unique_ptr<PmeTestHardwareContext>> PmeSolveTest::s_pmeTestHardwareContexts;
268
269 /*! \brief Test for PME solving */
270 TEST_P(PmeSolveTest, ReproducesOutputs)
271 {
272     EXPECT_NO_THROW_GMX(runTest());
273 }
274
275 /* Valid input instances */
276
277 //! A couple of valid inputs for boxes.
278 std::vector<Matrix3x3> const c_sampleBoxes{
279     // normal box
280     Matrix3x3{ { 8.0F, 0.0F, 0.0F, 0.0F, 3.4F, 0.0F, 0.0F, 0.0F, 2.0F } },
281     // triclinic box
282     Matrix3x3{ { 7.0F, 0.0F, 0.0F, 0.0F, 4.1F, 0.0F, 3.5F, 2.0F, 12.2F } },
283 };
284
285 //! A couple of valid inputs for grid sizes
286 std::vector<IVec> const c_sampleGridSizes{ IVec{ 16, 12, 28 }, IVec{ 9, 7, 23 } };
287
288 //! Moved out from instantiations for readability
289 const auto c_inputBoxes = ::testing::ValuesIn(c_sampleBoxes);
290 //! Moved out from instantiations for readability
291 const auto c_inputGridSizes = ::testing::ValuesIn(c_sampleGridSizes);
292
293 //! 2 sample complex grids - only non-zero values have to be listed
294 std::vector<SparseComplexGridValuesInput> const c_sampleGrids{
295     SparseComplexGridValuesInput{
296             { IVec{ 0, 0, 0 }, t_complex{ 3.5F, 6.7F } },
297             { IVec{ 7, 0, 0 }, t_complex{ -2.5F, -0.7F } },
298             { IVec{ 3, 5, 7 }, t_complex{ -0.006F, 1e-8F } },
299             { IVec{ 3, 1, 2 }, t_complex{ 0.6F, 7.9F } },
300             { IVec{ 6, 2, 4 }, t_complex{ 30.1F, 2.45F } },
301     },
302     SparseComplexGridValuesInput{
303             { IVec{ 0, 4, 0 }, t_complex{ 0.0F, 0.3F } },
304             { IVec{ 4, 2, 7 }, t_complex{ 13.76F, -40.0F } },
305             { IVec{ 0, 6, 7 }, t_complex{ 3.6F, 0.0F } },
306             { IVec{ 2, 5, 10 }, t_complex{ 3.6F, 10.65F } },
307     }
308 };
309
310 //! Moved out from instantiations for readability
311 const auto c_inputGrids = ::testing::ValuesIn(c_sampleGrids);
312 //! Moved out from instantiations for readability
313 const auto c_inputEpsilon_r = ::testing::Values(1.2);
314 //! Moved out from instantiations for readability
315 const auto c_inputEwaldCoeff_q = ::testing::Values(2.0);
316 //! Moved out from instantiations for readability
317 const auto c_inputEwaldCoeff_lj = ::testing::Values(0.7);
318 //! Moved out from instantiations for readability
319 const auto c_inputMethods = ::testing::Values(PmeSolveAlgorithm::Coulomb, PmeSolveAlgorithm::LennardJones);
320
321 //! Instantiation of the PME solving test
322 INSTANTIATE_TEST_SUITE_P(SaneInput,
323                          PmeSolveTest,
324                          ::testing::Combine(c_inputBoxes,
325                                             c_inputGridSizes,
326                                             c_inputGrids,
327                                             c_inputEpsilon_r,
328                                             c_inputEwaldCoeff_q,
329                                             c_inputEwaldCoeff_lj,
330                                             c_inputMethods));
331
332 //! A few more instances to check that different ewaldCoeff_q actually affects results of the Coulomb solver
333 INSTANTIATE_TEST_SUITE_P(DifferentEwaldCoeffQ,
334                          PmeSolveTest,
335                          ::testing::Combine(c_inputBoxes,
336                                             c_inputGridSizes,
337                                             c_inputGrids,
338                                             c_inputEpsilon_r,
339                                             ::testing::Values(0.4),
340                                             c_inputEwaldCoeff_lj,
341                                             ::testing::Values(PmeSolveAlgorithm::Coulomb)));
342
343 //! A few more instances to check that different ewaldCoeff_lj actually affects results of the Lennard-Jones solver.
344 //! The value has to be approximately larger than 1 / (box dimensions) to have a meaningful output grid.
345 //! Previous value of 0.3 caused one of the grid cells to be less or greater than GMX_FLOAT_MIN, depending on the architecture.
346 INSTANTIATE_TEST_SUITE_P(DifferentEwaldCoeffLJ,
347                          PmeSolveTest,
348                          ::testing::Combine(c_inputBoxes,
349                                             c_inputGridSizes,
350                                             c_inputGrids,
351                                             c_inputEpsilon_r,
352                                             c_inputEwaldCoeff_q,
353                                             ::testing::Values(2.35),
354                                             ::testing::Values(PmeSolveAlgorithm::LennardJones)));
355
356 //! A few more instances to check that different epsilon_r actually affects results of all solvers
357 INSTANTIATE_TEST_SUITE_P(DifferentEpsilonR,
358                          PmeSolveTest,
359                          ::testing::Combine(c_inputBoxes,
360                                             c_inputGridSizes,
361                                             c_inputGrids,
362                                             testing::Values(1.9),
363                                             c_inputEwaldCoeff_q,
364                                             c_inputEwaldCoeff_lj,
365                                             c_inputMethods));
366
367 } // namespace
368 } // namespace test
369 } // namespace gmx