Move t_forcetable and table enums to forcetable header
[alexxy/gromacs.git] / src / gromacs / tables / forcetable.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2014,2015,2016,2017 by the GROMACS development team.
5  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 #ifndef GMX_TABLES_FORCETABLE_H
37 #define GMX_TABLES_FORCETABLE_H
38
39 /*! \libinternal \file
40  * \brief
41  * Old routines for table generation (will eventually be replaced)
42  *
43  * \inlibraryapi
44  * \author Erik Lindahl <erik.lindahl@gmail.com>
45  * \ingroup module_tables
46  */
47
48 #include <cstdio>
49
50 #include <memory>
51 #include <vector>
52
53 #include "gromacs/utility/alignedallocator.h"
54 #include "gromacs/utility/real.h"
55
56 struct EwaldCorrectionTables;
57 struct bondedtable_t;
58 struct interaction_const_t;
59
60 /*! \brief Flag to select user tables for make_tables */
61 #define GMX_MAKETABLES_FORCEUSER (1 << 0)
62 /*! \brief Flag to only make 1,4 pair tables for make_tables */
63 #define GMX_MAKETABLES_14ONLY (1 << 1)
64
65 //! \brief The types of interactions contained in the table
66 enum class TableInteraction : int
67 {
68     VdwRepulsionVdwDispersion,
69     ElectrostaticVdwRepulsionVdwDispersion,
70     Count
71 };
72
73 /* Different formats for table data. Cubic spline tables are typically stored
74  * with the four Y,F,G,H intermediate values (check tables.c for format), which
75  * makes it easy to load with a single 4-way SIMD instruction too.
76  * Linear tables only need one value per table point, or two if both V and F
77  * are calculated. However, with SIMD instructions this makes the loads unaligned,
78  * and in that case we store the data as F, D=F(i+1)-F(i), V, and then a blank value,
79  * which again makes it possible to load as a single instruction.
80  */
81 enum class TableFormat : int
82 {
83     CubicsplineYfgh,
84     Count
85 };
86
87 //! \internal \brief Structure describing the data in a single table
88 struct t_forcetable
89 {
90     t_forcetable(TableInteraction interaction, TableFormat format);
91
92     ~t_forcetable();
93
94     //! Types of interactions stored in this table
95     TableInteraction interaction;
96     //! Interpolation type and data format
97     TableFormat format;
98     //! range of the table
99     real interactionRange;
100     //! n+1 is the number of table points
101     int numTablePoints;
102     //! distance (nm) between two table points
103     real scale;
104     //! The actual table data
105     std::vector<real, gmx::AlignedAllocator<real>> data;
106
107     /* Some information about the table layout. This can also be derived from the interpolation
108      * type and the table interactions, but it is convenient to have here for sanity checks, and it
109      * makes it much easier to access the tables in the nonbonded kernels when we can set the data
110      * from variables. It is always true that stride = formatsize*ninteractions
111      */
112
113     //! Number of fp variables for each table point (1 for F, 2 for VF, 4 for YFGH, etc.), only YFGH is implemented
114     static constexpr int formatsize = 4;
115     //! Number of interactions in table, 1 for coul-only, 3 for coul+rep+disp.
116     int numInteractions;
117     //! Distance to next table point (number of fp variables per table point in total)
118     int stride;
119 };
120
121 /*! \brief Enumerated type to describe the interaction types in a table */
122 enum
123 {
124     etiCOUL, //!< Coulomb
125     etiLJ6,  //!< Dispersion
126     etiLJ12, //!< Repulsion
127     etiNR    //!< Total number of interaction types
128 };
129
130 /*! \brief Function pointer to calculate the grid contribution for coulomb/LJ
131  *
132  * Used to tell table_spline3_fill_ewald_lr whether it
133  * should calculate the grid contribution for electrostatics or LJ.
134  */
135 typedef double (*real_space_grid_contribution_computer)(double, double);
136
137
138 /*! \brief Construct tables with the Ewald long-range force interaction
139  *
140  * Creates and fills tables of numPoints points with the spacing
141  * set to 1/tableScaling with the Ewald long-range (mesh) force.
142  * There are three separate tables with format F, V, FDV0.
143  * This function interpolates the Ewald mesh potential contribution
144  * with coefficient beta using a quadratic spline.
145  * The force is then be interpolated linearly.
146  *
147  * \param numPoints     Number of points in the tables
148  * \param tableScaling  1/spacing, units 1/nm
149  * \param beta          Ewald splitting parameter, units 1/nm
150  * \param v_lr          Pointer to function calculating real-space grid contribution
151  * \returns a set of Ewald correction tables
152  */
153 EwaldCorrectionTables generateEwaldCorrectionTables(int    numPoints,
154                                                     double tableScaling,
155                                                     real   beta,
156                                                     real_space_grid_contribution_computer v_lr);
157
158 /*! \brief Compute scaling for the Ewald quadratic spline tables.
159  *
160  * \param ic                     Pointer to interaction constant structure
161  * \param generateCoulombTables  Take the spacing for Coulomb Ewald corrections into account
162  * \param generateVdwTables      Take the spacing for Van der Waals Ewald corrections into account
163  * \return The scaling factor in units 1/nm
164  */
165 real ewald_spline3_table_scale(const interaction_const_t& ic, bool generateCoulombTables, bool generateVdwTables);
166
167 /*! \brief Return the real space grid contribution for Ewald
168  *
169  *  \param beta  Ewald splitting parameter
170  *  \param r     Distance for which to calculate the real-space contrib
171  *  \return      Real space grid contribution for Ewald electrostatics
172  */
173 double v_q_ewald_lr(double beta, double r);
174
175 /*! \brief Return the real space grid contribution for LJ-Ewald
176  *
177  *  \param beta  Ewald splitting parameter
178  *  \param r     Distance for which to calculate the real-space contrib
179  *  \return      Real space grid contribution for Ewald Lennard-Jones interaction
180  */
181 double v_lj_ewald_lr(double beta, double r);
182
183 /*! \brief Return tables for inner loops.
184  *
185  * \param fp     Log file pointer
186  * \param ic     Non-bonded interaction constants
187  * \param fn     File name from which to read user tables
188  * \param rtab   Largest interaction distance to tabulate
189  * \param flags  Flags to select table settings
190  *
191  * \return Pointer to inner loop table structure
192  */
193 std::unique_ptr<t_forcetable>
194 make_tables(FILE* fp, const interaction_const_t* ic, const char* fn, real rtab, int flags);
195
196 /*! \brief Return a table for bonded interactions,
197  *
198  * \param  fplog   Pointer to log file
199  * \param  fn      File name
200  * \param  angle   Type of angle: bonds 0, angles 1, dihedrals 2
201  * \return New bonded table datatype
202  */
203 bondedtable_t make_bonded_table(FILE* fplog, const char* fn, int angle);
204
205 /*! \brief Construct and return tabulated dispersion and repulsion interactions
206  *
207  * This table can be used to compute long-range dispersion corrections.
208  * Returns pointer owning nothing when tabfn=nullptr.
209  */
210 std::unique_ptr<t_forcetable>
211 makeDispersionCorrectionTable(FILE* fp, const interaction_const_t* ic, real rtab, const char* tabfn);
212
213 #endif /* GMX_TABLES_FORCETABLE_H */