Remove gmx custom fixed int (e.g. gmx_int64_t) types
[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, 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/utility/arrayref.h"
55 #include "gromacs/utility/unique_cptr.h"
56
57 #include "testhardwarecontexts.h"
58
59 namespace gmx
60 {
61 namespace test
62 {
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>using SparseGridValuesInput = std::map<IVec, ValueType>;
85 //! Non-zero real grid values
86 typedef SparseGridValuesInput<real> SparseRealGridValuesInput;
87 //! Non-zero complex grid values
88 typedef SparseGridValuesInput<t_complex> SparseComplexGridValuesInput;
89 //! Non-zero grid values for test output; keys are string representations of the cells' 3d indices (IVec); this allows for better sorting.
90 template<typename ValueType>using SparseGridValuesOutput = std::map<std::string, ValueType>;
91 //! Non-zero real grid values
92 typedef SparseGridValuesOutput<real> SparseRealGridValuesOutput;
93 //! Non-zero complex grid values
94 typedef SparseGridValuesOutput<t_complex> SparseComplexGridValuesOutput;
95 //! TODO: make proper C++ matrix for the whole Gromacs, get rid of this
96 typedef std::array<real, DIM * DIM> Matrix3x3;
97 //! PME solver type
98 enum class PmeSolveAlgorithm
99 {
100     Coulomb,
101     LennardJones,
102 };
103 //! PME solver results - reciprocal energy and virial
104 typedef std::tuple<real, Matrix3x3> PmeSolveOutput;
105
106 // Misc.
107
108 //! Tells if this generally valid PME input is supported for this mode
109 bool pmeSupportsInputForMode(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 //! Simple PME initialization (no atom data)
121 PmeSafePointer pmeInitEmpty(const t_inputrec *inputRec,
122                             CodePath mode = CodePath::CPU,
123                             gmx_device_info_t *gpuInfo = nullptr,
124                             PmeGpuProgramHandle pmeGpuProgram = nullptr,
125                             const Matrix3x3 &box = {{1.0f, 0.0f, 0.0f, 0.0f, 1.0f, 0.0f, 0.0f, 0.0f, 1.0f}},
126                             real ewaldCoeff_q = 0.0f, real ewaldCoeff_lj = 0.0f);
127 //! PME initialization with atom data and system box
128 PmeSafePointer pmeInitAtoms(const t_inputrec         *inputRec,
129                             CodePath                  mode,
130                             gmx_device_info_t        *gpuInfo,
131                             PmeGpuProgramHandle       pmeGpuProgram,
132                             const CoordinatesVector  &coordinates,
133                             const ChargesVector      &charges,
134                             const Matrix3x3          &box
135                             );
136 //! PME spline computation and charge spreading
137 void pmePerformSplineAndSpread(gmx_pme_t *pme, CodePath mode,
138                                bool computeSplines, bool spreadCharges);
139 //! PME solving
140 void pmePerformSolve(const gmx_pme_t *pme, CodePath mode,
141                      PmeSolveAlgorithm method, real cellVolume,
142                      GridOrdering gridOrdering, bool computeEnergyAndVirial);
143 //! PME force gathering
144 void pmePerformGather(gmx_pme_t *pme, CodePath mode,
145                       PmeForceOutputHandling inputTreatment, ForcesVector &forces);
146 //! PME test finalization before fetching the outputs
147 void pmeFinalizeTest(const gmx_pme_t *pme, CodePath mode);
148
149 // PME state setters
150
151 //! Setting atom spline values or derivatives to be used in spread/gather
152 void pmeSetSplineData(const gmx_pme_t *pme, CodePath mode,
153                       const SplineParamsDimVector &splineValues, PmeSplineDataType type, int dimIndex);
154 //! Setting gridline indices be used in spread/gather
155 void pmeSetGridLineIndices(const gmx_pme_t *pme, CodePath mode,
156                            const GridLineIndicesVector &gridLineIndices);
157 //! Setting real grid to be used in gather
158 void pmeSetRealGrid(const gmx_pme_t *pme, CodePath mode,
159                     const SparseRealGridValuesInput &gridValues);
160 void pmeSetComplexGrid(const gmx_pme_t *pme, CodePath mode, GridOrdering gridOrdering,
161                        const SparseComplexGridValuesInput &gridValues);
162
163 // PME state getters
164
165 //! Getting the single dimension's spline values or derivatives
166 SplineParamsDimVector pmeGetSplineData(const gmx_pme_t *pme, CodePath mode,
167                                        PmeSplineDataType type, int dimIndex);
168 //! Getting the gridline indices
169 GridLineIndicesVector pmeGetGridlineIndices(const gmx_pme_t *pme, CodePath mode);
170 //! Getting the real grid (spreading output of pmePerformSplineAndSpread())
171 SparseRealGridValuesOutput pmeGetRealGrid(const gmx_pme_t *pme, CodePath mode);
172 //! Getting the complex grid output of pmePerformSolve()
173 SparseComplexGridValuesOutput pmeGetComplexGrid(const gmx_pme_t *pme, CodePath mode,
174                                                 GridOrdering gridOrdering);
175 //! Getting the reciprocal energy and virial
176 PmeSolveOutput pmeGetReciprocalEnergyAndVirial(const gmx_pme_t *pme, CodePath mode,
177                                                PmeSolveAlgorithm method);
178 }
179 }
180
181 #endif