Use unique_ptr for tables
[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
52 #include "gromacs/utility/real.h"
53
54 struct EwaldCorrectionTables;
55 struct bondedtable_t;
56 struct interaction_const_t;
57 struct t_forcetable;
58
59 /*! \brief Flag to select user tables for make_tables */
60 #define GMX_MAKETABLES_FORCEUSER (1 << 0)
61 /*! \brief Flag to only make 1,4 pair tables for make_tables */
62 #define GMX_MAKETABLES_14ONLY (1 << 1)
63
64 /*! \brief Enumerated type to describe the interaction types in a table */
65 enum
66 {
67     etiCOUL, //!< Coulomb
68     etiLJ6,  //!< Dispersion
69     etiLJ12, //!< Repulsion
70     etiNR    //!< Total number of interaction types
71 };
72
73 /*! \brief Function pointer to calculate the grid contribution for coulomb/LJ
74  *
75  * Used to tell table_spline3_fill_ewald_lr whether it
76  * should calculate the grid contribution for electrostatics or LJ.
77  */
78 typedef double (*real_space_grid_contribution_computer)(double, double);
79
80
81 /*! \brief Construct tables with the Ewald long-range force interaction
82  *
83  * Creates and fills tables of numPoints points with the spacing
84  * set to 1/tableScaling with the Ewald long-range (mesh) force.
85  * There are three separate tables with format F, V, FDV0.
86  * This function interpolates the Ewald mesh potential contribution
87  * with coefficient beta using a quadratic spline.
88  * The force is then be interpolated linearly.
89  *
90  * \param numPoints     Number of points in the tables
91  * \param tableScaling  1/spacing, units 1/nm
92  * \param beta          Ewald splitting parameter, units 1/nm
93  * \param v_lr          Pointer to function calculating real-space grid contribution
94  * \returns a set of Ewald correction tables
95  */
96 EwaldCorrectionTables generateEwaldCorrectionTables(int    numPoints,
97                                                     double tableScaling,
98                                                     real   beta,
99                                                     real_space_grid_contribution_computer v_lr);
100
101 /*! \brief Compute scaling for the Ewald quadratic spline tables.
102  *
103  * \param ic                     Pointer to interaction constant structure
104  * \param generateCoulombTables  Take the spacing for Coulomb Ewald corrections into account
105  * \param generateVdwTables      Take the spacing for Van der Waals Ewald corrections into account
106  * \return The scaling factor in units 1/nm
107  */
108 real ewald_spline3_table_scale(const interaction_const_t& ic, bool generateCoulombTables, bool generateVdwTables);
109
110 /*! \brief Return the real space grid contribution for Ewald
111  *
112  *  \param beta  Ewald splitting parameter
113  *  \param r     Distance for which to calculate the real-space contrib
114  *  \return      Real space grid contribution for Ewald electrostatics
115  */
116 double v_q_ewald_lr(double beta, double r);
117
118 /*! \brief Return the real space grid contribution for LJ-Ewald
119  *
120  *  \param beta  Ewald splitting parameter
121  *  \param r     Distance for which to calculate the real-space contrib
122  *  \return      Real space grid contribution for Ewald Lennard-Jones interaction
123  */
124 double v_lj_ewald_lr(double beta, double r);
125
126 /*! \brief Return tables for inner loops.
127  *
128  * \param fp     Log file pointer
129  * \param ic     Non-bonded interaction constants
130  * \param fn     File name from which to read user tables
131  * \param rtab   Largest interaction distance to tabulate
132  * \param flags  Flags to select table settings
133  *
134  * \return Pointer to inner loop table structure
135  */
136 std::unique_ptr<t_forcetable>
137 make_tables(FILE* fp, const interaction_const_t* ic, const char* fn, real rtab, int flags);
138
139 /*! \brief Return a table for bonded interactions,
140  *
141  * \param  fplog   Pointer to log file
142  * \param  fn      File name
143  * \param  angle   Type of angle: bonds 0, angles 1, dihedrals 2
144  * \return New bonded table datatype
145  */
146 bondedtable_t make_bonded_table(FILE* fplog, const char* fn, int angle);
147
148 /*! \brief Construct and return tabulated dispersion and repulsion interactions
149  *
150  * This table can be used to compute long-range dispersion corrections.
151  * Returns pointer owning nothing when tabfn=nullptr.
152  */
153 std::unique_ptr<t_forcetable>
154 makeDispersionCorrectionTable(FILE* fp, const interaction_const_t* ic, real rtab, const char* tabfn);
155
156 #endif /* GMX_TABLES_FORCETABLE_H */