Move testing code to test directory
[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/arrayref.h"
54 #include "gromacs/utility/unique_cptr.h"
55
56 namespace gmx
57 {
58 namespace test
59 {
60
61 // Forward declaration
62 enum class CodePath;
63
64 // Convenience typedefs
65 //! A safe pointer type for PME.
66 typedef gmx::unique_cptr<gmx_pme_t, gmx_pme_destroy> PmeSafePointer;
67 //! Charges
68 typedef ArrayRef<const real> ChargesVector;
69 //! Coordinates
70 typedef std::vector<RVec> CoordinatesVector;
71 //! Forces
72 typedef ArrayRef<RVec> ForcesVector;
73 //! Gridline indices
74 typedef ArrayRef<const IVec> GridLineIndicesVector;
75 /*! \brief Spline parameters (theta or dtheta).
76  * A reference to a single dimension's spline data; this means (atomCount * pmeOrder) values or derivatives.
77  */
78 typedef ArrayRef<const real> SplineParamsDimVector;
79 /*! \brief Spline parameters (theta or dtheta) in all 3 dimensions
80  */
81 typedef std::array<SplineParamsDimVector, DIM> SplineParamsVector;
82
83 //! Non-zero grid values for test input; keys are 3d indices (IVec)
84 template<typename ValueType>
85 using SparseGridValuesInput = std::map<IVec, ValueType>;
86 //! Non-zero real grid values
87 typedef SparseGridValuesInput<real> SparseRealGridValuesInput;
88 //! Non-zero complex grid values
89 typedef SparseGridValuesInput<t_complex> SparseComplexGridValuesInput;
90 //! Non-zero grid values for test output; keys are string representations of the cells' 3d indices (IVec); this allows for better sorting.
91 template<typename ValueType>
92 using SparseGridValuesOutput = std::map<std::string, ValueType>;
93 //! Non-zero real grid values
94 typedef SparseGridValuesOutput<real> SparseRealGridValuesOutput;
95 //! Non-zero complex grid values
96 typedef SparseGridValuesOutput<t_complex> SparseComplexGridValuesOutput;
97 //! TODO: make proper C++ matrix for the whole Gromacs, get rid of this
98 typedef std::array<real, DIM * DIM> Matrix3x3;
99 //! PME solver type
100 enum class PmeSolveAlgorithm
101 {
102     Coulomb,
103     LennardJones,
104 };
105
106 // Misc.
107
108 //! Tells if this generally valid PME input is supported for this mode
109 bool pmeSupportsInputForMode(const gmx_hw_info_t& hwinfo, const t_inputrec* inputRec, CodePath mode);
110
111 //! Spline moduli are computed in double precision, so they're very good in single precision
112 constexpr int64_t c_splineModuliSinglePrecisionUlps = 1;
113 /*! \brief For double precision checks, the recursive interpolation
114  * and use of trig functions in make_dft_mod require a lot more flops,
115  * and thus opportunity for deviation between implementations. */
116 uint64_t getSplineModuliDoublePrecisionUlps(int splineOrder);
117
118 // PME stages
119
120 //! PME initialization
121 PmeSafePointer pmeInitWrapper(const t_inputrec*        inputRec,
122                               CodePath                 mode,
123                               const DeviceInformation* deviceInfo,
124                               const PmeGpuProgram*     pmeGpuProgram,
125                               const Matrix3x3&         box,
126                               real                     ewaldCoeff_q  = 1.0F,
127                               real                     ewaldCoeff_lj = 1.0F);
128 //! Simple PME initialization (no atom data)
129 PmeSafePointer pmeInitEmpty(const t_inputrec*        inputRec,
130                             CodePath                 mode,
131                             const DeviceInformation* deviceInfo,
132                             const PmeGpuProgram*     pmeGpuProgram,
133                             const Matrix3x3&         box,
134                             real                     ewaldCoeff_q,
135                             real                     ewaldCoeff_lj);
136 //! Simple PME initialization based on inputrec only
137 PmeSafePointer pmeInitEmpty(const t_inputrec* inputRec);
138 //! Make a GPU state-propagator manager
139 std::unique_ptr<StatePropagatorDataGpu> makeStatePropagatorDataGpu(const gmx_pme_t&     pme,
140                                                                    const DeviceContext& deviceContext);
141 //! PME initialization with atom data and system box
142 void pmeInitAtoms(gmx_pme_t*               pme,
143                   StatePropagatorDataGpu*  stateGpu,
144                   CodePath                 mode,
145                   const CoordinatesVector& coordinates,
146                   const ChargesVector&     charges);
147 //! PME spline computation and charge spreading
148 void pmePerformSplineAndSpread(gmx_pme_t* pme, CodePath mode, bool computeSplines, bool spreadCharges);
149 //! PME solving
150 void pmePerformSolve(const gmx_pme_t*  pme,
151                      CodePath          mode,
152                      PmeSolveAlgorithm method,
153                      real              cellVolume,
154                      GridOrdering      gridOrdering,
155                      bool              computeEnergyAndVirial);
156 //! PME force gathering
157 void pmePerformGather(gmx_pme_t*    pme,
158                       CodePath      mode,
159                       ForcesVector& forces); //NOLINT(google-runtime-references)
160 //! PME test finalization before fetching the outputs
161 void pmeFinalizeTest(const gmx_pme_t* pme, CodePath mode);
162
163 // PME state setters
164
165 //! Setting atom spline values or derivatives to be used in spread/gather
166 void pmeSetSplineData(const gmx_pme_t*             pme,
167                       CodePath                     mode,
168                       const SplineParamsDimVector& splineValues,
169                       PmeSplineDataType            type,
170                       int                          dimIndex);
171
172 //! Setting gridline indices be used in spread/gather
173 void pmeSetGridLineIndices(gmx_pme_t* pme, CodePath mode, const GridLineIndicesVector& gridLineIndices);
174 //! Setting real grid to be used in gather
175 void pmeSetRealGrid(const gmx_pme_t* pme, CodePath mode, const SparseRealGridValuesInput& gridValues);
176 void pmeSetComplexGrid(const gmx_pme_t*                    pme,
177                        CodePath                            mode,
178                        GridOrdering                        gridOrdering,
179                        const SparseComplexGridValuesInput& gridValues);
180
181 // PME state getters
182
183 //! Getting the single dimension's spline values or derivatives
184 SplineParamsDimVector pmeGetSplineData(const gmx_pme_t* pme, CodePath mode, PmeSplineDataType type, int dimIndex);
185 //! Getting the gridline indices
186 GridLineIndicesVector pmeGetGridlineIndices(const gmx_pme_t* pme, CodePath mode);
187 //! Getting the real grid (spreading output of pmePerformSplineAndSpread())
188 SparseRealGridValuesOutput pmeGetRealGrid(const gmx_pme_t* pme, CodePath mode);
189 //! Getting the complex grid output of pmePerformSolve()
190 SparseComplexGridValuesOutput pmeGetComplexGrid(const gmx_pme_t* pme, CodePath mode, GridOrdering gridOrdering);
191 //! Getting the reciprocal energy and virial
192 PmeOutput pmeGetReciprocalEnergyAndVirial(const gmx_pme_t* pme, CodePath mode, PmeSolveAlgorithm method);
193 } // namespace test
194 } // namespace gmx
195
196 #endif