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