Merge branch release-2016
[alexxy/gromacs.git] / src / gromacs / ewald / tests / pmesplinespreadtest.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 spline computation and charge spreading 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
64 //! PME spline and spread code path being tested
65 enum class PmeSplineAndSpreadOptions
66 {
67     SplineOnly,
68     SpreadOnly,
69     SplineAndSpreadUnified
70 };
71
72 /*! \brief Convenience typedef of input parameters - unit cell box, PME interpolation order, grid dimensions,
73  * particle coordinates, particle charges
74  * TODO: consider inclusion of local grid offsets/sizes or PME nodes counts to test the PME DD
75  */
76 typedef std::tuple<Matrix3x3, int, IVec, CoordinatesVector, ChargesVector> SplineAndSpreadInputParameters;
77
78 /*! \brief Test fixture for testing both atom spline parameter computation and charge spreading.
79  * These 2 stages of PME are tightly coupled in the code.
80  */
81 class PmeSplineAndSpreadTest : public ::testing::TestWithParam<SplineAndSpreadInputParameters>
82 {
83     public:
84         //! Default constructor
85         PmeSplineAndSpreadTest() = default;
86         //! The test
87         void runTest()
88         {
89             /* Getting the input */
90             Matrix3x3               box;
91             int                     pmeOrder;
92             IVec                    gridSize;
93             CoordinatesVector       coordinates;
94             ChargesVector           charges;
95
96             std::tie(box, pmeOrder, gridSize, coordinates, charges) = GetParam();
97             const size_t atomCount = coordinates.size();
98
99             /* Storing the input where it's needed */
100             t_inputrec inputRec;
101             inputRec.nkx         = gridSize[XX];
102             inputRec.nky         = gridSize[YY];
103             inputRec.nkz         = gridSize[ZZ];
104             inputRec.pme_order   = pmeOrder;
105             inputRec.coulombtype = eelPME;
106
107             TestReferenceData                                      refData;
108
109             const std::map<CodePath, std::string>                  modesToTest   = {{CodePath::CPU, "CPU"}};
110             const std::map<PmeSplineAndSpreadOptions, std::string> optionsToTest = {{PmeSplineAndSpreadOptions::SplineAndSpreadUnified, "spline computation and charge spreading (fused)"},
111                                                                                     {PmeSplineAndSpreadOptions::SplineOnly, "spline computation"},
112                                                                                     {PmeSplineAndSpreadOptions::SpreadOnly, "charge spreading"}};
113             for (const auto &mode : modesToTest)
114             {
115                 for (const auto &option : optionsToTest)
116                 {
117                     /* Describing the test uniquely in case it fails */
118
119                     SCOPED_TRACE(formatString("Testing %s with %s for PME grid size %d %d %d"
120                                               ", order %d, %zu atoms",
121                                               option.second.c_str(), mode.second.c_str(),
122                                               gridSize[XX], gridSize[YY], gridSize[ZZ],
123                                               pmeOrder,
124                                               atomCount));
125
126                     /* Running the test */
127
128                     PmeSafePointer pmeSafe = pmeInitWithAtoms(&inputRec, coordinates, charges, box);
129
130                     const bool     computeSplines = (option.first == PmeSplineAndSpreadOptions::SplineOnly) || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
131                     const bool     spreadCharges  = (option.first == PmeSplineAndSpreadOptions::SpreadOnly) || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
132
133                     if (!computeSplines)
134                     {
135                         // Here we should set up the results of the spline computation so that the spread can run.
136                         // What is lazy and works is running the separate spline so that it will set it up for us:
137                         pmePerformSplineAndSpread(pmeSafe.get(), mode.first, true, false);
138                         // We know that it is tested in another iteration.
139                         // TODO: Clean alternative: read and set the reference gridline indices, spline params
140                     }
141
142                     pmePerformSplineAndSpread(pmeSafe.get(), mode.first, computeSplines, spreadCharges);
143
144                     /* Outputs correctness check */
145                     /* All tolerances were picked empirically for single precision on CPU */
146
147                     TestReferenceChecker rootChecker(refData.rootChecker());
148
149                     const auto           maxGridSize              = std::max(std::max(gridSize[XX], gridSize[YY]), gridSize[ZZ]);
150                     const auto           ulpToleranceSplineValues = 2 * (pmeOrder - 2) * maxGridSize;
151                     /* 2 is empiric, the rest follows from the amount of operations */
152
153                     if (computeSplines)
154                     {
155                         const char *dimString[] = { "X", "Y", "Z" };
156
157                         /* Spline values */
158                         SCOPED_TRACE(formatString("Testing spline values with tolerance of %ld", ulpToleranceSplineValues));
159                         TestReferenceChecker splineValuesChecker(rootChecker.checkCompound("Splines", "Values"));
160                         splineValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceSplineValues));
161                         for (int i = 0; i < DIM; i++)
162                         {
163                             auto splineValuesDim = pmeGetSplineData(pmeSafe.get(), mode.first, PmeSplineDataType::Values, i);
164                             splineValuesChecker.checkSequence(splineValuesDim.begin(), splineValuesDim.end(), dimString[i]);
165                         }
166
167                         /* Spline derivatives */
168                         const auto ulpToleranceSplineDerivatives = 4 * ulpToleranceSplineValues;
169                         /* 4 is just a wild guess since the derivatives are deltas of neighbor spline values which could differ greatly */
170                         SCOPED_TRACE(formatString("Testing spline derivatives with tolerance of %ld", ulpToleranceSplineDerivatives));
171                         TestReferenceChecker splineDerivativesChecker(rootChecker.checkCompound("Splines", "Derivatives"));
172                         splineDerivativesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceSplineDerivatives));
173                         for (int i = 0; i < DIM; i++)
174                         {
175                             auto splineDerivativesDim = pmeGetSplineData(pmeSafe.get(), mode.first, PmeSplineDataType::Derivatives, i);
176                             splineDerivativesChecker.checkSequence(splineDerivativesDim.begin(), splineDerivativesDim.end(), dimString[i]);
177                         }
178
179                         /* Particle gridline indices */
180                         auto gridLineIndices = pmeGetGridlineIndices(pmeSafe.get(), mode.first);
181                         rootChecker.checkSequence(gridLineIndices.begin(), gridLineIndices.end(), "Gridline indices");
182                     }
183
184                     if (spreadCharges)
185                     {
186                         /* The wrapped grid */
187                         SparseRealGridValuesOutput nonZeroGridValues = pmeGetRealGrid(pmeSafe.get(), mode.first);
188                         TestReferenceChecker       gridValuesChecker(rootChecker.checkCompound("NonZeroGridValues", "RealSpaceGrid"));
189                         const auto                 ulpToleranceGrid = 2 * ulpToleranceSplineValues * (int)(ceil(sqrt(atomCount)));
190                         /* 2 is empiric; sqrt(atomCount) assumes all the input charges may spread onto the same cell */
191                         SCOPED_TRACE(formatString("Testing grid values with tolerance of %ld", ulpToleranceGrid));
192                         gridValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceGrid));
193                         for (const auto &point : nonZeroGridValues)
194                         {
195                             gridValuesChecker.checkReal(point.second, point.first.c_str());
196                         }
197                     }
198                 }
199             }
200         }
201 };
202
203
204 /*! \brief Test for PME B-spline moduli computation */
205 TEST_P(PmeSplineAndSpreadTest, ReproducesOutputs)
206 {
207     EXPECT_NO_THROW(runTest());
208 }
209
210 /* Valid input instances */
211
212 //! A couple of valid inputs for boxes.
213 static std::vector<Matrix3x3> const sampleBoxes
214 {
215     // normal box
216     Matrix3x3 {{
217                    8.0f, 0.0f, 0.0f,
218                    0.0f, 3.4f, 0.0f,
219                    0.0f, 0.0f, 2.0f
220                }},
221     // triclinic box
222     Matrix3x3 {{
223                    7.0f, 0.0f, 0.0f,
224                    0.0f, 4.1f, 0.0f,
225                    3.5f, 2.0f, 12.2f
226                }},
227 };
228
229 //! A couple of valid inputs for grid sizes.
230 static std::vector<IVec> const sampleGridSizes
231 {
232     IVec {
233         16, 12, 14
234     },
235     IVec {
236         19, 17, 11
237     }
238 };
239
240 //! Random charges
241 static std::vector<real> const sampleChargesFull
242 {
243     4.95f, 3.11f, 3.97f, 1.08f, 2.09f, 1.1f, 4.13f, 3.31f, 2.8f, 5.83f, 5.09f, 6.1f, 2.86f, 0.24f, 5.76f, 5.19f, 0.72f
244 };
245 //! 1 charge
246 static auto const sampleCharges1 = ChargesVector::fromVector(sampleChargesFull.begin(), sampleChargesFull.begin() + 1);
247 //! 2 charges
248 static auto const sampleCharges2 = ChargesVector::fromVector(sampleChargesFull.begin() + 1, sampleChargesFull.begin() + 3);
249 //! 13 charges
250 static auto const sampleCharges13 = ChargesVector::fromVector(sampleChargesFull.begin() + 3, sampleChargesFull.begin() + 16);
251
252 //! Random coordinate vectors
253 static CoordinatesVector const sampleCoordinatesFull
254 {
255     {
256         5.59f, 1.37f, 0.95f
257     }, {
258         16.0f, 1.02f, 0.22f      // 2 box lengths in x
259     }, {
260         0.034f, 1.65f, 0.22f
261     }, {
262         0.33f, 0.92f, 1.56f
263     }, {
264         1.16f, 0.75f, 0.39f
265     }, {
266         0.5f, 1.63f, 1.14f
267     }, {
268         16.0001f, 1.52f, 1.19f  // > 2 box lengths in x
269     }, {
270         1.43f, 1.1f, 4.1f       // > 2 box lengths in z
271     }, {
272         -1.08f, 1.19f, 0.08f    // negative x
273     }, {
274         1.6f, 0.93f, 0.53f
275     }, {
276         1.32f, -1.48f, 0.16f    // negative y
277     }, {
278         0.87f, 0.0f, 0.33f
279     }, {
280         0.95f, 7.7f, -0.48f     // > 2 box lengths in y, negative z
281     }, {
282         1.23f, 0.91f, 0.68f
283     }, {
284         0.19f, 1.45f, 0.94f
285     }, {
286         1.28f, 0.46f, 0.38f
287     }, {
288         1.21f, 0.23f, 1.0f
289     }
290 };
291 //! 1 coordinate vector
292 static CoordinatesVector const sampleCoordinates1(sampleCoordinatesFull.begin(), sampleCoordinatesFull.begin() + 1);
293 //! 2 coordinate vectors
294 static CoordinatesVector const sampleCoordinates2(sampleCoordinatesFull.begin() + 1, sampleCoordinatesFull.begin() + 3);
295 //! 13 coordinate vectors
296 static CoordinatesVector const sampleCoordinates13(sampleCoordinatesFull.begin() + 3, sampleCoordinatesFull.begin() + 16);
297
298 //! moved out from instantiantions for readability
299 auto inputBoxes     = ::testing::ValuesIn(sampleBoxes);
300 //! moved out from instantiantions for readability
301 auto inputPmeOrders = ::testing::Range(3, 5 + 1);
302 //! moved out from instantiantions for readability
303 auto inputGridSizes = ::testing::ValuesIn(sampleGridSizes);
304
305 /*! \brief Instantiation of the PME spline computation test with valid input and 1 atom */
306 INSTANTIATE_TEST_CASE_P(SaneInput1, PmeSplineAndSpreadTest, ::testing::Combine(inputBoxes, inputPmeOrders, inputGridSizes,
307                                                                                    ::testing::Values(sampleCoordinates1),
308                                                                                    ::testing::Values(sampleCharges1)
309                                                                                ));
310 /*! \brief Instantiation of the PME spline computation test with valid input and 2 atoms */
311 INSTANTIATE_TEST_CASE_P(SaneInput2, PmeSplineAndSpreadTest, ::testing::Combine(inputBoxes, inputPmeOrders, inputGridSizes,
312                                                                                    ::testing::Values(sampleCoordinates2),
313                                                                                    ::testing::Values(sampleCharges2)
314                                                                                ));
315 /*! \brief Instantiation of the PME spline computation test with valid input and 13 atoms */
316 INSTANTIATE_TEST_CASE_P(SaneInput13, PmeSplineAndSpreadTest, ::testing::Combine(inputBoxes, inputPmeOrders, inputGridSizes,
317                                                                                     ::testing::Values(sampleCoordinates13),
318                                                                                     ::testing::Values(sampleCharges13)
319                                                                                 ));
320 }
321 }
322 }