Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / tables / forcetable.h
blobf98327beafdb69e98976e800eb7397a5c59ae6da
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2014,2015,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 #ifndef GMX_TABLES_FORCETABLE_H
36 #define GMX_TABLES_FORCETABLE_H
38 /*! \libinternal \file
39 * \brief
40 * Old routines for table generation (will eventually be replaced)
42 * \inlibraryapi
43 * \author Erik Lindahl <erik.lindahl@gmail.com>
44 * \ingroup module_tables
47 #include <cstdio>
49 #include <memory>
51 #include "gromacs/mdtypes/fcdata.h"
52 #include "gromacs/mdtypes/forcerec.h"
53 #include "gromacs/mdtypes/interaction_const.h"
54 #include "gromacs/utility/real.h"
56 /*! \brief Flag to select user tables for make_tables */
57 #define GMX_MAKETABLES_FORCEUSER (1<<0)
58 /*! \brief Flag to only make 1,4 pair tables for make_tables */
59 #define GMX_MAKETABLES_14ONLY (1<<1)
61 /*! \brief Enumerated type to describe the interaction types in a table */
62 enum {
63 etiCOUL, //!< Coulomb
64 etiLJ6, //!< Dispersion
65 etiLJ12, //!< Repulsion
66 etiNR //!< Total number of interaction types
69 /*! \brief Function pointer to calculate the grid contribution for coulomb/LJ
71 * Used to tell table_spline3_fill_ewald_lr whether it
72 * should calculate the grid contribution for electrostatics or LJ.
74 typedef double (*real_space_grid_contribution_computer)(double, double);
77 /*! \brief Construct tables with the Ewald long-range force interaction
79 * Creates and fills tables of numPoints points with the spacing
80 * set to 1/tableScaling with the Ewald long-range (mesh) force.
81 * There are three separate tables with format F, V, FDV0.
82 * This function interpolates the Ewald mesh potential contribution
83 * with coefficient beta using a quadratic spline.
84 * The force is then be interpolated linearly.
86 * \param numPoints Number of points in the tables
87 * \param tableScaling 1/spacing, units 1/nm
88 * \param beta Ewald splitting parameter, units 1/nm
89 * \param v_lr Pointer to function calculating real-space grid contribution
90 * \returns a set of Ewald correction tables
92 EwaldCorrectionTables
93 generateEwaldCorrectionTables(int numPoints,
94 double tableScaling,
95 real beta,
96 real_space_grid_contribution_computer v_lr);
98 /*! \brief Compute scaling for the Ewald quadratic spline tables.
100 * \param ic Pointer to interaction constant structure
101 * \param generateCoulombTables Take the spacing for Coulomb Ewald corrections into account
102 * \param generateVdwTables Take the spacing for Van der Waals Ewald corrections into account
103 * \return The scaling factor in units 1/nm
105 real ewald_spline3_table_scale(const interaction_const_t &ic,
106 bool generateCoulombTables,
107 bool generateVdwTables);
109 /*! \brief Return the real space grid contribution for Ewald
111 * \param beta Ewald splitting parameter
112 * \param r Distance for which to calculate the real-space contrib
113 * \return Real space grid contribution for Ewald electrostatics
115 double v_q_ewald_lr(double beta, double r);
117 /*! \brief Return the real space grid contribution for LJ-Ewald
119 * \param beta Ewald splitting parameter
120 * \param r Distance for which to calculate the real-space contrib
121 * \return Real space grid contribution for Ewald Lennard-Jones interaction
123 double v_lj_ewald_lr(double beta, double r);
125 /*! \brief Return tables for inner loops.
127 * \param fp Log file pointer
128 * \param ic Non-bonded interaction constants
129 * \param fn File name from which to read user tables
130 * \param rtab Largest interaction distance to tabulate
131 * \param flags Flags to select table settings
133 * \return Pointer to inner loop table structure
135 t_forcetable *make_tables(FILE *fp,
136 const interaction_const_t *ic,
137 const char *fn, real rtab, int flags);
139 /*! \brief Return a table for bonded interactions,
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
146 bondedtable_t make_bonded_table(FILE *fplog, const char *fn, int angle);
148 /*! \brief Construct and return tabulated dispersion and repulsion interactions
150 * This table can be used to compute long-range dispersion corrections.
151 * Returns pointer owning nothing when tabfn=nullptr.
153 std::unique_ptr<t_forcetable>
154 makeDispersionCorrectionTable(FILE *fp, const interaction_const_t *ic,
155 real rtab, const char *tabfn);
157 #endif /* GMX_TABLES_FORCETABLE_H */