2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017, 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 <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/utility/arrayref.h"
55 #include "gromacs/utility/unique_cptr.h"
57 #include "testhardwarecontexts.h"
64 // Convenience typedefs
65 //! A safe pointer type for PME.
66 typedef gmx::unique_cptr<gmx_pme_t, gmx_pme_destroy> PmeSafePointer;
68 typedef ArrayRef<const real> ChargesVector;
70 typedef std::vector<RVec> CoordinatesVector;
72 typedef ArrayRef<RVec> ForcesVector;
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.
78 typedef ArrayRef<const real> SplineParamsDimVector;
79 /*! \brief Spline parameters (theta or dtheta) in all 3 dimensions
81 typedef std::array<SplineParamsDimVector, DIM> SplineParamsVector;
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 gathering input forces treatment
98 enum class PmeGatherInputHandling
104 enum class PmeSolveAlgorithm
109 //! PME solver results - reciprocal energy and virial
110 typedef std::tuple<real, Matrix3x3> PmeSolveOutput;
114 //! Tells if this generally valid PME input is supported for this mode
115 bool pmeSupportsInputForMode(const t_inputrec *inputRec, CodePath mode);
117 //! Returns tolerance of anything directly influenced by B-Splines (relaxed for double precision)
118 class FloatingPointTolerance getSplineTolerance(gmx_int64_t toleranceUlps);
122 // TODO: currently PME initializations do not store CodePath. They probably should (unless we would need mixed CPU-GPU execution?).
123 //! Simple PME initialization (no atom data)
124 PmeSafePointer pmeInitEmpty(const t_inputrec *inputRec,
125 CodePath mode = CodePath::CPU,
126 gmx_device_info_t *gpuInfo = nullptr,
127 const Matrix3x3 &box = {{1.0f, 0.0f, 0.0f, 0.0f, 1.0f, 0.0f, 0.0f, 0.0f, 1.0f}},
128 real ewaldCoeff_q = 0.0f, real ewaldCoeff_lj = 0.0f);
129 //! PME initialization with atom data and system box
130 PmeSafePointer pmeInitAtoms(const t_inputrec *inputRec,
132 gmx_device_info_t *gpuInfo,
133 const CoordinatesVector &coordinates,
134 const ChargesVector &charges,
137 //! PME spline computation and charge spreading
138 void pmePerformSplineAndSpread(gmx_pme_t *pme, CodePath mode,
139 bool computeSplines, bool spreadCharges);
141 void pmePerformSolve(const gmx_pme_t *pme, CodePath mode,
142 PmeSolveAlgorithm method, real cellVolume,
143 GridOrdering gridOrdering, bool computeEnergyAndVirial);
144 //! PME force gathering
145 void pmePerformGather(gmx_pme_t *pme, CodePath mode,
146 PmeGatherInputHandling inputTreatment, ForcesVector &forces);
147 //! PME test finalization before fetching the outputs
148 void pmeFinalizeTest(const gmx_pme_t *pme, CodePath mode);
152 //! Setting atom spline values or derivatives to be used in spread/gather
153 void pmeSetSplineData(const gmx_pme_t *pme, CodePath mode,
154 const SplineParamsDimVector &splineValues, PmeSplineDataType type, int dimIndex);
155 //! Setting gridline indices be used in spread/gather
156 void pmeSetGridLineIndices(const gmx_pme_t *pme, CodePath mode,
157 const GridLineIndicesVector &gridLineIndices);
158 //! Setting real grid to be used in gather
159 void pmeSetRealGrid(const gmx_pme_t *pme, CodePath mode,
160 const SparseRealGridValuesInput &gridValues);
161 void pmeSetComplexGrid(const gmx_pme_t *pme, CodePath mode, GridOrdering gridOrdering,
162 const SparseComplexGridValuesInput &gridValues);
166 //! Getting the single dimension's spline values or derivatives
167 SplineParamsDimVector pmeGetSplineData(const gmx_pme_t *pme, CodePath mode,
168 PmeSplineDataType type, int dimIndex);
169 //! Getting the gridline indices
170 GridLineIndicesVector pmeGetGridlineIndices(const gmx_pme_t *pme, CodePath mode);
171 //! Getting the real grid (spreading output of pmePerformSplineAndSpread())
172 SparseRealGridValuesOutput pmeGetRealGrid(const gmx_pme_t *pme, CodePath mode);
173 //! Getting the complex grid output of pmePerformSolve()
174 SparseComplexGridValuesOutput pmeGetComplexGrid(const gmx_pme_t *pme, CodePath mode,
175 GridOrdering gridOrdering);
176 //! Getting the reciprocal energy and virial
177 PmeSolveOutput pmeGetReciprocalEnergyAndVirial(const gmx_pme_t *pme, CodePath mode,
178 PmeSolveAlgorithm method);