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 "gromacs/ewald/pme.h"
50 #include "gromacs/math/gmxcomplex.h"
51 #include "gromacs/utility/arrayref.h"
52 #include "gromacs/utility/unique_cptr.h"
61 // Convenience typedefs
62 //! A safe pointer type for PME.
63 typedef gmx::unique_cptr<gmx_pme_t, gmx_pme_destroy> PmeSafePointer;
65 typedef ConstArrayRef<real> ChargesVector;
67 typedef std::vector<RVec> CoordinatesVector;
69 typedef ArrayRef<RVec> ForcesVector;
71 typedef ConstArrayRef<IVec> GridLineIndicesVector;
72 //! Type of spline data
73 enum class PmeSplineDataType
76 Derivatives, // dtheta
78 /*! \brief Spline parameters (theta or dtheta).
79 * A reference to a single dimension's spline data; this means (atomCount * pmeOrder) values or derivatives.
81 typedef ConstArrayRef<real> SplineParamsDimVector;
82 /*! \brief Spline parameters (theta or dtheta) in all 3 dimensions
84 typedef std::array<SplineParamsDimVector, DIM> SplineParamsVector;
86 //! Non-zero grid values for test input; keys are 3d indices (IVec)
87 template<typename ValueType>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>using SparseGridValuesOutput = std::map<std::string, ValueType>;
94 //! Non-zero real grid values
95 typedef SparseGridValuesOutput<real> SparseRealGridValuesOutput;
96 //! Non-zero complex grid values
97 typedef SparseGridValuesOutput<t_complex> SparseComplexGridValuesOutput;
98 //! TODO: make proper C++ matrix for the whole Gromacs, get rid of this
99 typedef std::array<real, DIM * DIM> Matrix3x3;
100 //! PME code path being tested
103 CPU, // serial CPU code
105 //! PME gathering input forces treatment
106 enum class PmeGatherInputHandling
112 enum class PmeSolveAlgorithm
117 //! PME grid dimension ordering (from major to minor)
118 enum class GridOrdering
123 //! PME solver results - reciprocal energy and virial
124 typedef std::tuple<real, Matrix3x3> PmeSolveOutput;
128 //! Simple PME initialization (no atom data)
129 PmeSafePointer pmeInitEmpty(const t_inputrec *inputRec,
130 const Matrix3x3 &box = {{1.0f, 0.0f, 0.0f, 0.0f, 1.0f, 0.0f, 0.0f, 0.0f, 1.0f}},
131 real ewaldCoeff_q = 0.0f, real ewaldCoeff_lj = 0.0f);
132 //! PME initialization with atom data and system box; lacks Ewald coefficients
133 PmeSafePointer pmeInitAtoms(const t_inputrec *inputRec,
134 const CoordinatesVector &coordinates,
135 const ChargesVector &charges,
138 //! PME spline computation and charge spreading
139 void pmePerformSplineAndSpread(gmx_pme_t *pme, CodePath mode,
140 bool computeSplines, bool spreadCharges);
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 PmeGatherInputHandling inputTreatment, ForcesVector &forces);
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);
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);