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