Apply re-formatting to C++ in src/ tree.
[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,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 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     PmeBSplineModuliTest() = default;
80     //! The whole logic being tested is contained here
81     void runTest()
82     {
83         /* Getting the input */
84         const BSplineModuliInputParameters parameters = GetParam();
85         int                                pmeOrder;
86         IVec                               gridSize;
87         ModuliType                         moduliType;
88         std::tie(gridSize, pmeOrder, moduliType) = parameters;
89
90         /* Describing the test in case it fails */
91         SCOPED_TRACE(formatString(
92                 "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],
96                 gridSize[YY],
97                 gridSize[ZZ]));
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.coulombtype = (moduliType == ModuliType::P3M) ? eelP3M_AD : eelPME;
105         inputRec.pme_order   = pmeOrder;
106
107         /* PME initialization call which checks the inputs and computes the B-spline moduli according to the grid sizes. */
108         PmeSafePointer pme = pmeInitEmpty(&inputRec);
109
110         /* Setting up the checker */
111         TestReferenceData    refData;
112         TestReferenceChecker checker(refData.rootChecker());
113         checker.setDefaultTolerance(relativeToleranceAsPrecisionDependentUlp(
114                 1.0, c_splineModuliSinglePrecisionUlps, getSplineModuliDoublePrecisionUlps(pmeOrder + 1)));
115
116         /* Perform a correctness check */
117         const char* dimString[] = { "X", "Y", "Z" };
118         for (int i = 0; i < DIM; i++)
119         {
120             checker.checkSequenceArray(gridSize[i], pme->bsp_mod[i], dimString[i]);
121         }
122     }
123 };
124
125 // Different aliases are needed to reuse the test class without instantiating tests more than once
126 //! Failure test alias
127 using PmeBSplineModuliFailureTest = PmeBSplineModuliTest;
128 TEST_P(PmeBSplineModuliFailureTest, Throws)
129 {
130     EXPECT_THROW_GMX(runTest(), InconsistentInputError);
131 }
132 //! Correctness test alias
133 using PmeBSplineModuliCorrectnessTest = PmeBSplineModuliTest;
134 TEST_P(PmeBSplineModuliCorrectnessTest, ReproducesValues)
135 {
136     EXPECT_NO_THROW_GMX(runTest());
137 }
138
139 /* Invalid input instances */
140
141 //! Sane interpolation order
142 const int sanePmeOrder = 4;
143 //! Sane grid size
144 const IVec saneGridSize = { 32, 25, 47 };
145 /*! \brief Hand-picked invalid input for the exception tests */
146 std::vector<BSplineModuliInputParameters> const invalidInputs{
147     /* Invalid grid sizes */
148     BSplineModuliInputParameters{ IVec{ -1, 10, 10 }, sanePmeOrder, ModuliType::P3M },
149     BSplineModuliInputParameters{ IVec{ 40, 0, 20 }, sanePmeOrder, ModuliType::P3M },
150     BSplineModuliInputParameters{ IVec{ 64, 2, 64 }, sanePmeOrder, ModuliType::PME },
151     /* Invalid interpolation orders */
152     BSplineModuliInputParameters{
153             saneGridSize,
154             8 + 1,
155             ModuliType::P3M // P3M only supports orders up to 8
156     },
157     BSplineModuliInputParameters{ saneGridSize, PME_ORDER_MAX + 1, ModuliType::PME },
158 };
159
160 /*! \brief Instantiation of the PME B-spline moduli creation test with invalid input */
161 INSTANTIATE_TEST_CASE_P(InsaneInput, PmeBSplineModuliFailureTest, ::testing::ValuesIn(invalidInputs));
162
163 /* Valid input instances */
164
165 //! A couple of valid inputs for grid sizes. It is good to test both even and odd dimensions.
166 std::vector<IVec> const sampleGridSizes{ IVec{ 64, 32, 64 }, IVec{ 57, 84, 29 } };
167
168 /*! \brief Instantiation of the PME B-spline moduli creation test with valid input - up to order of 8 */
169 INSTANTIATE_TEST_CASE_P(SaneInput1,
170                         PmeBSplineModuliCorrectnessTest,
171                         ::testing::Combine(::testing::ValuesIn(sampleGridSizes),
172                                            ::testing::Range(3, 8 + 1),
173                                            ::testing::Values(ModuliType::PME, ModuliType::P3M)));
174 } // namespace
175 } // namespace test
176 } // namespace gmx