clang-tidy: google tests applicable
[alexxy/gromacs.git] / src / gromacs / ewald / tests / pmebsplinetest.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 B-spline moduli computation 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/ewald/pme-internal.h"
50 #include "gromacs/math/vectypes.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/utility/stringutil.h"
53
54 #include "testutils/refdata.h"
55 #include "testutils/testasserts.h"
56
57 #include "pmetestcommon.h"
58
59 namespace gmx
60 {
61 namespace test
62 {
63 namespace
64 {
65 /*! \brief Moduli algorithm type */
66 enum class ModuliType
67 {
68     PME,
69     P3M
70 };
71
72 /*! \brief Convenience typedef of input parameters - grid dimensions, PME interpolation order, moduli type */
73 typedef std::tuple<IVec, int, ModuliType> BSplineModuliInputParameters;
74
75 /*! \brief Test fixture for testing PME B-spline moduli creation */
76 class PmeBSplineModuliTest : public ::testing::TestWithParam<BSplineModuliInputParameters>
77 {
78     public:
79         //! Default constructor
80         PmeBSplineModuliTest() = default;
81         //! The whole logic being tested is contained here
82         void runTest()
83         {
84             /* Getting the input */
85             const BSplineModuliInputParameters parameters = GetParam();
86             int         pmeOrder;
87             IVec        gridSize;
88             ModuliType  moduliType;
89             std::tie(gridSize, pmeOrder, moduliType) = parameters;
90
91             /* Describing the test in case it fails */
92             SCOPED_TRACE(formatString("Testing B-spline moduli creation (%s) for PME order %d, grid size %d %d %d",
93                                       (moduliType == ModuliType::P3M) ? "P3M" : "plain",
94                                       pmeOrder,
95                                       gridSize[XX], gridSize[YY], gridSize[ZZ]));
96
97             /* Storing the input where it's needed */
98             t_inputrec inputRec;
99             inputRec.nkx         = gridSize[XX];
100             inputRec.nky         = gridSize[YY];
101             inputRec.nkz         = gridSize[ZZ];
102             inputRec.coulombtype = (moduliType == ModuliType::P3M) ? eelP3M_AD : eelPME;
103             inputRec.pme_order   = pmeOrder;
104
105             /* PME initialization call which checks the inputs and computes the B-spline moduli according to the grid sizes. */
106             PmeSafePointer pme = pmeInitEmpty(&inputRec);
107
108             /* Setting up the checker */
109             TestReferenceData    refData;
110             TestReferenceChecker checker(refData.rootChecker());
111             checker.setDefaultTolerance(relativeToleranceAsPrecisionDependentUlp(1.0,
112                                                                                  c_splineModuliSinglePrecisionUlps,
113                                                                                  getSplineModuliDoublePrecisionUlps(pmeOrder+1)));
114
115             /* Perform a correctness check */
116             const char *dimString[] = { "X", "Y", "Z" };
117             for (int i = 0; i < DIM; i++)
118             {
119                 checker.checkSequenceArray(gridSize[i], pme->bsp_mod[i], dimString[i]);
120             }
121         }
122 };
123
124 // Different aliases are needed to reuse the test class without instantiating tests more than once
125 //! Failure test alias
126 using PmeBSplineModuliFailureTest = PmeBSplineModuliTest;
127 TEST_P(PmeBSplineModuliFailureTest, Throws)
128 {
129     EXPECT_THROW(runTest(), InconsistentInputError);
130 }
131 //! Correctness test alias
132 using PmeBSplineModuliCorrectnessTest = PmeBSplineModuliTest;
133 TEST_P(PmeBSplineModuliCorrectnessTest, ReproducesValues)
134 {
135     EXPECT_NO_THROW(runTest());
136 }
137
138 /* Invalid input instances */
139
140 //! Sane interpolation order
141 const int       sanePmeOrder = 4;
142 //! Sane grid size
143 const IVec      saneGridSize = {32, 25, 47};
144 /*! \brief Hand-picked invalid input for the exception tests */
145 std::vector<BSplineModuliInputParameters> const invalidInputs
146 {
147     /* Invalid grid sizes */
148     BSplineModuliInputParameters {
149         IVec {
150             -1, 10, 10
151         }, sanePmeOrder, ModuliType::P3M
152     },
153     BSplineModuliInputParameters {
154         IVec {
155             40, 0, 20
156         }, sanePmeOrder, ModuliType::P3M
157     },
158     BSplineModuliInputParameters {
159         IVec {
160             64, 2, 64
161         }, sanePmeOrder, ModuliType::PME
162     },
163     /* Invalid interpolation orders */
164     BSplineModuliInputParameters {
165         saneGridSize, 8 + 1, ModuliType::P3M  // P3M only supports orders up to 8
166     },
167     BSplineModuliInputParameters {
168         saneGridSize, PME_ORDER_MAX + 1, ModuliType::PME
169     },
170 };
171
172 /*! \brief Instantiation of the PME B-spline moduli creation test with invalid input */
173 INSTANTIATE_TEST_CASE_P(InsaneInput, PmeBSplineModuliFailureTest, ::testing::ValuesIn(invalidInputs));
174
175 /* Valid input instances */
176
177 //! A couple of valid inputs for grid sizes. It is good to test both even and odd dimensions.
178 std::vector<IVec> const sampleGridSizes
179 {
180     IVec {
181         64, 32, 64
182     },
183     IVec {
184         57, 84, 29
185     }
186 };
187
188 /*! \brief Instantiation of the PME B-spline moduli creation test with valid input - up to order of 8 */
189 INSTANTIATE_TEST_CASE_P(SaneInput1, PmeBSplineModuliCorrectnessTest, ::testing::Combine(
190                                     ::testing::ValuesIn(sampleGridSizes),
191                                     ::testing::Range(3, 8 + 1),
192                                     ::testing::Values(ModuliType::PME, ModuliType::P3M)
193                                 ));
194 }  // namespace
195 }  // namespace test
196 }  // namespace gmx