2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 * Describes common routines and types for PME tests.
39 * \author Aleksei Iupinov <a.yupinov@gmail.com>
40 * \ingroup module_ewald
42 #ifndef GMX_EWALD_PME_TEST_COMMON_H
43 #define GMX_EWALD_PME_TEST_COMMON_H
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"
58 class DeviceStreamManager;
65 // Forward declaration
68 // Convenience typedefs
69 //! A safe pointer type for PME.
70 typedef gmx::unique_cptr<gmx_pme_t, gmx_pme_destroy> PmeSafePointer;
72 typedef ArrayRef<const real> ChargesVector;
74 typedef std::vector<RVec> CoordinatesVector;
76 typedef ArrayRef<RVec> ForcesVector;
78 typedef ArrayRef<const IVec> GridLineIndicesVector;
79 /*! \brief Spline parameters (theta or dtheta).
80 * A reference to a single dimension's spline data; this means (atomCount * pmeOrder) values or derivatives.
82 typedef ArrayRef<const real> SplineParamsDimVector;
83 /*! \brief Spline parameters (theta or dtheta) in all 3 dimensions
85 typedef std::array<SplineParamsDimVector, DIM> SplineParamsVector;
87 //! Non-zero grid values for test input; keys are 3d indices (IVec)
88 template<typename ValueType>
89 using SparseGridValuesInput = std::map<IVec, ValueType>;
90 //! Non-zero real grid values
91 typedef SparseGridValuesInput<real> SparseRealGridValuesInput;
92 //! Non-zero complex grid values
93 typedef SparseGridValuesInput<t_complex> SparseComplexGridValuesInput;
94 //! Non-zero grid values for test output; keys are string representations of the cells' 3d indices (IVec); this allows for better sorting.
95 template<typename ValueType>
96 using SparseGridValuesOutput = std::map<std::string, ValueType>;
97 //! Non-zero real grid values
98 typedef SparseGridValuesOutput<real> SparseRealGridValuesOutput;
99 //! Non-zero complex grid values
100 typedef SparseGridValuesOutput<t_complex> SparseComplexGridValuesOutput;
101 //! TODO: make proper C++ matrix for the whole Gromacs, get rid of this
102 typedef std::array<real, DIM * DIM> Matrix3x3;
104 enum class PmeSolveAlgorithm
112 //! Tells if this generally valid PME input is supported for this mode
113 bool pmeSupportsInputForMode(const gmx_hw_info_t& hwinfo, const t_inputrec* inputRec, CodePath mode);
115 //! Spline moduli are computed in double precision, so they're very good in single precision
116 constexpr int64_t c_splineModuliSinglePrecisionUlps = 1;
117 /*! \brief For double precision checks, the recursive interpolation
118 * and use of trig functions in make_dft_mod require a lot more flops,
119 * and thus opportunity for deviation between implementations. */
120 uint64_t getSplineModuliDoublePrecisionUlps(int splineOrder);
124 //! PME initialization
125 PmeSafePointer pmeInitWrapper(const t_inputrec* inputRec,
127 const DeviceContext* deviceContext,
128 const DeviceStream* deviceStream,
129 const PmeGpuProgram* pmeGpuProgram,
130 const Matrix3x3& box,
131 real ewaldCoeff_q = 1.0F,
132 real ewaldCoeff_lj = 1.0F);
133 //! Simple PME initialization (no atom data)
134 PmeSafePointer pmeInitEmpty(const t_inputrec* inputRec,
136 const DeviceContext* deviceContext,
137 const DeviceStream* deviceStream,
138 const PmeGpuProgram* pmeGpuProgram,
139 const Matrix3x3& box,
143 //! Simple PME initialization based on inputrec only
144 PmeSafePointer pmeInitEmpty(const t_inputrec* inputRec);
146 //! Make a GPU state-propagator manager
147 std::unique_ptr<StatePropagatorDataGpu> makeStatePropagatorDataGpu(const gmx_pme_t& pme,
148 const DeviceContext* deviceContext,
149 const DeviceStream* deviceStream);
150 //! PME initialization with atom data and system box
151 void pmeInitAtoms(gmx_pme_t* pme,
152 StatePropagatorDataGpu* stateGpu,
154 const CoordinatesVector& coordinates,
155 const ChargesVector& charges);
156 //! PME spline computation and charge spreading
157 void pmePerformSplineAndSpread(gmx_pme_t* pme, CodePath mode, bool computeSplines, bool spreadCharges);
159 void pmePerformSolve(const gmx_pme_t* pme,
161 PmeSolveAlgorithm method,
163 GridOrdering gridOrdering,
164 bool computeEnergyAndVirial);
165 //! PME force gathering
166 void pmePerformGather(gmx_pme_t* pme,
168 ForcesVector& forces); //NOLINT(google-runtime-references)
169 //! PME test finalization before fetching the outputs
170 void pmeFinalizeTest(const gmx_pme_t* pme, CodePath mode);
174 //! Setting atom spline values or derivatives to be used in spread/gather
175 void pmeSetSplineData(const gmx_pme_t* pme,
177 const SplineParamsDimVector& splineValues,
178 PmeSplineDataType type,
181 //! Setting gridline indices be used in spread/gather
182 void pmeSetGridLineIndices(gmx_pme_t* pme, CodePath mode, const GridLineIndicesVector& gridLineIndices);
183 //! Setting real grid to be used in gather
184 void pmeSetRealGrid(const gmx_pme_t* pme, CodePath mode, const SparseRealGridValuesInput& gridValues);
185 void pmeSetComplexGrid(const gmx_pme_t* pme,
187 GridOrdering gridOrdering,
188 const SparseComplexGridValuesInput& gridValues);
192 //! Getting the single dimension's spline values or derivatives
193 SplineParamsDimVector pmeGetSplineData(const gmx_pme_t* pme, CodePath mode, PmeSplineDataType type, int dimIndex);
194 //! Getting the gridline indices
195 GridLineIndicesVector pmeGetGridlineIndices(const gmx_pme_t* pme, CodePath mode);
196 //! Getting the real grid (spreading output of pmePerformSplineAndSpread())
197 SparseRealGridValuesOutput pmeGetRealGrid(const gmx_pme_t* pme, CodePath mode);
198 //! Getting the complex grid output of pmePerformSolve()
199 SparseComplexGridValuesOutput pmeGetComplexGrid(const gmx_pme_t* pme, CodePath mode, GridOrdering gridOrdering);
200 //! Getting the reciprocal energy and virial
201 PmeOutput pmeGetReciprocalEnergyAndVirial(const gmx_pme_t* pme, CodePath mode, PmeSolveAlgorithm method);