Refactor PME tests for better usability
[alexxy/gromacs.git] / src / gromacs / ewald / tests / pmetestcommon.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018,2019,2020,2021, 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  * Describes common routines and types for PME tests.
38  *
39  * \author Aleksei Iupinov <a.yupinov@gmail.com>
40  * \ingroup module_ewald
41  */
42 #ifndef GMX_EWALD_PME_TEST_COMMON_H
43 #define GMX_EWALD_PME_TEST_COMMON_H
44
45 #include <array>
46 #include <map>
47 #include <vector>
48
49 #include "gromacs/ewald/pme.h"
50 #include "gromacs/ewald/pme_gpu_internal.h"
51 #include "gromacs/math/gmxcomplex.h"
52 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
53 #include "gromacs/utility/unique_cptr.h"
54
55 #include "testutils/test_device.h"
56
57 namespace gmx
58 {
59
60 template<typename>
61 class ArrayRef;
62
63 namespace test
64 {
65
66 //! Hardware code path being tested
67 enum class CodePath : int
68 {
69     //! CPU code path
70     CPU,
71     //! GPU code path
72     GPU,
73     //! Total number of code paths
74     Count
75 };
76
77 // Convenience typedefs
78 //! A safe pointer type for PME.
79 typedef gmx::unique_cptr<gmx_pme_t, gmx_pme_destroy> PmeSafePointer;
80 //! Charges
81 typedef std::vector<real> ChargesVector;
82 //! Coordinates
83 typedef std::vector<RVec> CoordinatesVector;
84 //! Forces
85 typedef ArrayRef<RVec> ForcesVector;
86 //! Gridline indices
87 typedef std::vector<IVec> GridLineIndicesVector;
88 /*! \brief Spline parameters (theta or dtheta).
89  * A reference to a single dimension's spline data; this means (atomCount * pmeOrder) values or derivatives.
90  */
91 typedef ArrayRef<const real> SplineParamsDimVector;
92 /*! \brief Spline parameters (theta or dtheta) in all 3 dimensions
93  */
94 typedef std::array<SplineParamsDimVector, DIM> SplineParamsVector;
95
96 //! Non-zero grid values for test input; keys are 3d indices (IVec)
97 template<typename ValueType>
98 using SparseGridValuesInput = std::map<IVec, ValueType>;
99 //! Non-zero real grid values
100 typedef SparseGridValuesInput<real> SparseRealGridValuesInput;
101 //! Non-zero complex grid values
102 typedef SparseGridValuesInput<t_complex> SparseComplexGridValuesInput;
103 //! Non-zero grid values for test output; keys are string representations of the cells' 3d indices (IVec); this allows for better sorting.
104 template<typename ValueType>
105 using SparseGridValuesOutput = std::map<std::string, ValueType>;
106 //! Non-zero real grid values
107 typedef SparseGridValuesOutput<real> SparseRealGridValuesOutput;
108 //! Non-zero complex grid values
109 typedef SparseGridValuesOutput<t_complex> SparseComplexGridValuesOutput;
110 //! TODO: make proper C++ matrix for the whole Gromacs, get rid of this
111 typedef std::array<real, DIM * DIM> Matrix3x3;
112 //! PME solver type
113 enum class PmeSolveAlgorithm : int
114 {
115     //! Coulomb electrostatics
116     Coulomb,
117     //! Lennard-Jones
118     LennardJones,
119     //! Total number of solvers
120     Count
121 };
122
123 // Misc.
124
125 //! Tells if this generally valid PME input is supported for this mode
126 bool pmeSupportsInputForMode(const gmx_hw_info_t& hwinfo, const t_inputrec* inputRec, CodePath mode);
127
128 //! Spline moduli are computed in double precision, so they're very good in single precision
129 constexpr int64_t c_splineModuliSinglePrecisionUlps = 1;
130 /*! \brief For double precision checks, the recursive interpolation
131  * and use of trig functions in make_dft_mod require a lot more flops,
132  * and thus opportunity for deviation between implementations. */
133 uint64_t getSplineModuliDoublePrecisionUlps(int splineOrder);
134
135 // PME stages
136
137 //! PME initialization
138 PmeSafePointer pmeInitWrapper(const t_inputrec*    inputRec,
139                               CodePath             mode,
140                               const DeviceContext* deviceContext,
141                               const DeviceStream*  deviceStream,
142                               const PmeGpuProgram* pmeGpuProgram,
143                               const Matrix3x3&     box,
144                               real                 ewaldCoeff_q  = 1.0F,
145                               real                 ewaldCoeff_lj = 1.0F);
146
147 //! Simple PME initialization based on inputrec only
148 PmeSafePointer pmeInitEmpty(const t_inputrec* inputRec);
149
150 //! Make a GPU state-propagator manager
151 std::unique_ptr<StatePropagatorDataGpu> makeStatePropagatorDataGpu(const gmx_pme_t&     pme,
152                                                                    const DeviceContext* deviceContext,
153                                                                    const DeviceStream* deviceStream);
154 //! PME initialization with atom data and system box
155 void pmeInitAtoms(gmx_pme_t*               pme,
156                   StatePropagatorDataGpu*  stateGpu,
157                   CodePath                 mode,
158                   const CoordinatesVector& coordinates,
159                   const ChargesVector&     charges);
160 //! PME spline computation and charge spreading
161 void pmePerformSplineAndSpread(gmx_pme_t* pme, CodePath mode, bool computeSplines, bool spreadCharges);
162 //! PME solving
163 void pmePerformSolve(const gmx_pme_t*  pme,
164                      CodePath          mode,
165                      PmeSolveAlgorithm method,
166                      real              cellVolume,
167                      GridOrdering      gridOrdering,
168                      bool              computeEnergyAndVirial);
169 //! PME force gathering
170 void pmePerformGather(gmx_pme_t*    pme,
171                       CodePath      mode,
172                       ForcesVector& forces); //NOLINT(google-runtime-references)
173 //! PME test finalization before fetching the outputs
174 void pmeFinalizeTest(const gmx_pme_t* pme, CodePath mode);
175
176 // PME state setters
177
178 //! Setting atom spline values or derivatives to be used in spread/gather
179 void pmeSetSplineData(const gmx_pme_t*             pme,
180                       CodePath                     mode,
181                       const SplineParamsDimVector& splineValues,
182                       PmeSplineDataType            type,
183                       int                          dimIndex);
184
185 //! Setting gridline indices be used in spread/gather
186 void pmeSetGridLineIndices(gmx_pme_t* pme, CodePath mode, const GridLineIndicesVector& gridLineIndices);
187 //! Setting real grid to be used in gather
188 void pmeSetRealGrid(const gmx_pme_t* pme, CodePath mode, const SparseRealGridValuesInput& gridValues);
189 void pmeSetComplexGrid(const gmx_pme_t*                    pme,
190                        CodePath                            mode,
191                        GridOrdering                        gridOrdering,
192                        const SparseComplexGridValuesInput& gridValues);
193
194 // PME state getters
195
196 //! Getting the single dimension's spline values or derivatives
197 SplineParamsDimVector pmeGetSplineData(const gmx_pme_t* pme, CodePath mode, PmeSplineDataType type, int dimIndex);
198 //! Getting the gridline indices
199 GridLineIndicesVector pmeGetGridlineIndices(const gmx_pme_t* pme, CodePath mode);
200 //! Getting the real grid (spreading output of pmePerformSplineAndSpread())
201 SparseRealGridValuesOutput pmeGetRealGrid(const gmx_pme_t* pme, CodePath mode);
202 //! Getting the complex grid output of pmePerformSolve()
203 SparseComplexGridValuesOutput pmeGetComplexGrid(const gmx_pme_t* pme, CodePath mode, GridOrdering gridOrdering);
204 //! Getting the reciprocal energy and virial
205 PmeOutput pmeGetReciprocalEnergyAndVirial(const gmx_pme_t* pme, CodePath mode, PmeSolveAlgorithm method);
206
207 struct PmeTestHardwareContext
208 {
209     //! Hardware path for the code being tested.
210     CodePath codePath_;
211     //! Returns a human-readable context description line
212     std::string description() const;
213     //! Pointer to the global test hardware device (if on GPU)
214     TestDevice* testDevice_ = nullptr;
215     //! PME GPU program if needed
216     PmeGpuProgramStorage pmeGpuProgram_ = nullptr;
217     // Constructor for CPU context
218     PmeTestHardwareContext();
219     // Constructor for GPU context
220     explicit PmeTestHardwareContext(TestDevice* testDevice);
221
222     //! Get the code path
223     CodePath codePath() const { return codePath_; }
224     //! Get the PME GPU program
225     const PmeGpuProgram* pmeGpuProgram() const
226     {
227         return codePath() == CodePath::GPU ? pmeGpuProgram_.get() : nullptr;
228     }
229
230     const DeviceContext* deviceContext() const
231     {
232         return codePath() == CodePath::GPU ? &testDevice_->deviceContext() : nullptr;
233     }
234
235     const DeviceStream* deviceStream() const
236     {
237         return codePath() == CodePath::GPU ? &testDevice_->deviceStream() : nullptr;
238     }
239
240     //! Activate the context (set the device)
241     void activate() const;
242 };
243
244 std::vector<std::unique_ptr<PmeTestHardwareContext>> createPmeTestHardwareContextList();
245
246 //! A couple of valid inputs for boxes.
247 extern const std::map<std::string, Matrix3x3> c_inputBoxes;
248
249 //! Valid PME orders for testing
250 extern std::vector<int> c_inputPmeOrders;
251
252 } // namespace test
253 } // namespace gmx
254
255 #endif // GMX_EWALD_PME_TEST_COMMON_H