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