StatePropagatorDataGpu object to manage GPU forces, positions and velocities buffers
[gromacs.git] / src / gromacs / ewald / tests / pmetestcommon.h
blobb51551c29a6041c558c74435010682209d96ef77
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
35 /*! \internal \file
36 * \brief
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
45 #include <array>
46 #include <map>
47 #include <vector>
49 #include <gtest/gtest.h>
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"
58 #include "testhardwarecontexts.h"
60 namespace gmx
62 namespace test
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.
79 typedef ArrayRef<const real> SplineParamsDimVector;
80 /*! \brief Spline parameters (theta or dtheta) in all 3 dimensions
82 typedef std::array<SplineParamsDimVector, DIM> SplineParamsVector;
84 //! Non-zero grid values for test input; keys are 3d indices (IVec)
85 template<typename ValueType>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>using SparseGridValuesOutput = std::map<std::string, ValueType>;
92 //! Non-zero real grid values
93 typedef SparseGridValuesOutput<real> SparseRealGridValuesOutput;
94 //! Non-zero complex grid values
95 typedef SparseGridValuesOutput<t_complex> SparseComplexGridValuesOutput;
96 //! TODO: make proper C++ matrix for the whole Gromacs, get rid of this
97 typedef std::array<real, DIM * DIM> Matrix3x3;
98 //! PME solver type
99 enum class PmeSolveAlgorithm
101 Coulomb,
102 LennardJones,
105 // Misc.
107 //! Tells if this generally valid PME input is supported for this mode
108 bool pmeSupportsInputForMode(const gmx_hw_info_t &hwinfo,
109 const t_inputrec *inputRec,
110 CodePath mode);
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);
119 // PME stages
121 //! Simple PME initialization (no atom data)
122 PmeSafePointer pmeInitEmpty(const t_inputrec *inputRec,
123 CodePath mode = CodePath::CPU,
124 const gmx_device_info_t *gpuInfo = nullptr,
125 PmeGpuProgramHandle pmeGpuProgram = nullptr,
126 const Matrix3x3 &box = {{1.0F, 0.0F, 0.0F, 0.0F, 1.0F, 0.0F, 0.0F, 0.0F, 1.0F}},
127 real ewaldCoeff_q = 0.0F, real ewaldCoeff_lj = 0.0F);
128 //! PME initialization with atom data and system box
129 PmeSafePointer pmeInitAtoms(const t_inputrec *inputRec,
130 CodePath mode,
131 const gmx_device_info_t *gpuInfo,
132 PmeGpuProgramHandle pmeGpuProgram,
133 const CoordinatesVector &coordinates,
134 const ChargesVector &charges,
135 const Matrix3x3 &box,
136 std::shared_ptr<StatePropagatorDataGpu> stateGpu
138 //! PME spline computation and charge spreading
139 void pmePerformSplineAndSpread(gmx_pme_t *pme, CodePath mode,
140 bool computeSplines, bool spreadCharges);
141 //! PME solving
142 void pmePerformSolve(const gmx_pme_t *pme, CodePath mode,
143 PmeSolveAlgorithm method, real cellVolume,
144 GridOrdering gridOrdering, bool computeEnergyAndVirial);
145 //! PME force gathering
146 void pmePerformGather(gmx_pme_t *pme, CodePath mode,
147 PmeForceOutputHandling inputTreatment, ForcesVector &forces); //NOLINT(google-runtime-references)
148 //! PME test finalization before fetching the outputs
149 void pmeFinalizeTest(const gmx_pme_t *pme, CodePath mode);
151 // PME state setters
153 //! Setting atom spline values or derivatives to be used in spread/gather
154 void pmeSetSplineData(const gmx_pme_t *pme, CodePath mode,
155 const SplineParamsDimVector &splineValues, PmeSplineDataType type, int dimIndex);
156 //! Setting gridline indices be used in spread/gather
157 void pmeSetGridLineIndices(gmx_pme_t *pme, CodePath mode,
158 const GridLineIndicesVector &gridLineIndices);
159 //! Setting real grid to be used in gather
160 void pmeSetRealGrid(const gmx_pme_t *pme, CodePath mode,
161 const SparseRealGridValuesInput &gridValues);
162 void pmeSetComplexGrid(const gmx_pme_t *pme, CodePath mode, GridOrdering gridOrdering,
163 const SparseComplexGridValuesInput &gridValues);
165 // PME state getters
167 //! Getting the single dimension's spline values or derivatives
168 SplineParamsDimVector pmeGetSplineData(const gmx_pme_t *pme, CodePath mode,
169 PmeSplineDataType type, int dimIndex);
170 //! Getting the gridline indices
171 GridLineIndicesVector pmeGetGridlineIndices(const gmx_pme_t *pme, CodePath mode);
172 //! Getting the real grid (spreading output of pmePerformSplineAndSpread())
173 SparseRealGridValuesOutput pmeGetRealGrid(const gmx_pme_t *pme, CodePath mode);
174 //! Getting the complex grid output of pmePerformSolve()
175 SparseComplexGridValuesOutput pmeGetComplexGrid(const gmx_pme_t *pme, CodePath mode,
176 GridOrdering gridOrdering);
177 //! Getting the reciprocal energy and virial
178 PmeOutput pmeGetReciprocalEnergyAndVirial(const gmx_pme_t *pme, CodePath mode,
179 PmeSolveAlgorithm method);
180 } // namespace test
181 } // namespace gmx
183 #endif