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