Use common device context header for OpenCL
[gromacs.git] / src / gromacs / ewald / tests / pmebsplinetest.cpp
blobb6a4e9e35f81b69e24f16180b26b36bfcabcd898
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
35 /*! \internal \file
36 * \brief
37 * Implements PME B-spline moduli computation tests.
39 * \author Aleksei Iupinov <a.yupinov@gmail.com>
40 * \ingroup module_ewald
43 #include "gmxpre.h"
45 #include <string>
47 #include <gmock/gmock.h>
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"
54 #include "testutils/refdata.h"
55 #include "testutils/testasserts.h"
57 #include "pmetestcommon.h"
59 namespace gmx
61 namespace test
63 namespace
65 /*! \brief Moduli algorithm type */
66 enum class ModuliType
68 PME,
69 P3M
72 /*! \brief Convenience typedef of input parameters - grid dimensions, PME interpolation order, moduli type */
73 typedef std::tuple<IVec, int, ModuliType> BSplineModuliInputParameters;
75 /*! \brief Test fixture for testing PME B-spline moduli creation */
76 class PmeBSplineModuliTest : public ::testing::TestWithParam<BSplineModuliInputParameters>
78 public:
79 PmeBSplineModuliTest() = default;
80 //! The whole logic being tested is contained here
81 void runTest()
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;
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", pmeOrder, gridSize[XX],
94 gridSize[YY], gridSize[ZZ]));
96 /* Storing the input where it's needed */
97 t_inputrec inputRec;
98 inputRec.nkx = gridSize[XX];
99 inputRec.nky = gridSize[YY];
100 inputRec.nkz = gridSize[ZZ];
101 inputRec.coulombtype = (moduliType == ModuliType::P3M) ? eelP3M_AD : eelPME;
102 inputRec.pme_order = pmeOrder;
104 /* PME initialization call which checks the inputs and computes the B-spline moduli according to the grid sizes. */
105 PmeSafePointer pme = pmeInitEmpty(&inputRec);
107 /* Setting up the checker */
108 TestReferenceData refData;
109 TestReferenceChecker checker(refData.rootChecker());
110 checker.setDefaultTolerance(relativeToleranceAsPrecisionDependentUlp(
111 1.0, c_splineModuliSinglePrecisionUlps, getSplineModuliDoublePrecisionUlps(pmeOrder + 1)));
113 /* Perform a correctness check */
114 const char* dimString[] = { "X", "Y", "Z" };
115 for (int i = 0; i < DIM; i++)
117 checker.checkSequenceArray(gridSize[i], pme->bsp_mod[i], dimString[i]);
122 // Different aliases are needed to reuse the test class without instantiating tests more than once
123 //! Failure test alias
124 using PmeBSplineModuliFailureTest = PmeBSplineModuliTest;
125 TEST_P(PmeBSplineModuliFailureTest, Throws)
127 EXPECT_THROW_GMX(runTest(), InconsistentInputError);
129 //! Correctness test alias
130 using PmeBSplineModuliCorrectnessTest = PmeBSplineModuliTest;
131 TEST_P(PmeBSplineModuliCorrectnessTest, ReproducesValues)
133 EXPECT_NO_THROW_GMX(runTest());
136 /* Invalid input instances */
138 //! Sane interpolation order
139 const int sanePmeOrder = 4;
140 //! Sane grid size
141 const IVec saneGridSize = { 32, 25, 47 };
142 /*! \brief Hand-picked invalid input for the exception tests */
143 std::vector<BSplineModuliInputParameters> const invalidInputs{
144 /* Invalid grid sizes */
145 BSplineModuliInputParameters{ IVec{ -1, 10, 10 }, sanePmeOrder, ModuliType::P3M },
146 BSplineModuliInputParameters{ IVec{ 40, 0, 20 }, sanePmeOrder, ModuliType::P3M },
147 BSplineModuliInputParameters{ IVec{ 64, 2, 64 }, sanePmeOrder, ModuliType::PME },
148 /* Invalid interpolation orders */
149 BSplineModuliInputParameters{
150 saneGridSize, 8 + 1, ModuliType::P3M // P3M only supports orders up to 8
152 BSplineModuliInputParameters{ saneGridSize, PME_ORDER_MAX + 1, ModuliType::PME },
155 /*! \brief Instantiation of the PME B-spline moduli creation test with invalid input */
156 INSTANTIATE_TEST_CASE_P(InsaneInput, PmeBSplineModuliFailureTest, ::testing::ValuesIn(invalidInputs));
158 /* Valid input instances */
160 //! A couple of valid inputs for grid sizes. It is good to test both even and odd dimensions.
161 std::vector<IVec> const sampleGridSizes{ IVec{ 64, 32, 64 }, IVec{ 57, 84, 29 } };
163 /*! \brief Instantiation of the PME B-spline moduli creation test with valid input - up to order of 8 */
164 INSTANTIATE_TEST_CASE_P(SaneInput1,
165 PmeBSplineModuliCorrectnessTest,
166 ::testing::Combine(::testing::ValuesIn(sampleGridSizes),
167 ::testing::Range(3, 8 + 1),
168 ::testing::Values(ModuliType::PME, ModuliType::P3M)));
169 } // namespace
170 } // namespace test
171 } // namespace gmx