2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, 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_MDLIB_DISPERSIONCORRECTION_H
36 #define GMX_MDLIB_DISPERSIONCORRECTION_H
43 #include "gromacs/math/vectypes.h"
46 struct interaction_const_t;
58 class DispersionCorrection
61 /*! \brief Constructor
63 * \param[in] mtop The global topology
64 * \param[in] inputrec The input record
65 * \param[in] useBuckingham True when Buckingham is used instead of LJ
66 * \param[in] numAtomTypes The number of non-bonded atom types
67 * \param[in] nonbondedForceParameters The LJ or Bham parameter matrix stored as a flat list
68 * \param[in] ic The nonbonded interaction parameters
69 * \param[in] tableFileName Table file name, should != nullptr (checked)
71 DispersionCorrection(const gmx_mtop_t& mtop,
72 const t_inputrec& inputrec,
75 gmx::ArrayRef<const real> nonbondedForceParameters,
76 const interaction_const_t& ic,
77 const char* tableFileName);
79 /*! \brief Print dispersion correction information to the log file
81 * \param[in] mdlog The MD logger
83 void print(const gmx::MDLogger& mdlog) const;
85 /*! \brief Computes and sets energy and virial correction parameters
87 * Sets all parameters that are affected by the cut-off and/or the
88 * LJ-Ewald coefficient. Should be called before calling calculate()
89 * and whenever interaction settings change, e.g. PME tuning.
91 * \param[in] ic The nonbonded interaction parameters
93 void setParameters(const interaction_const_t& ic);
96 * \brief Struct for returning all dispersion correction quantities
100 /*! \brief Correct the virial tensor for the missing dispersion
102 * \param[in,out] virialTensor The virial tensor to correct
104 void correctVirial(tensor virialTensor) const
106 for (int m = 0; m < DIM; m++)
108 virialTensor[m][m] += virial;
112 real virial = 0; //!< Scalar correction to the virial
113 real pressure = 0; //!< Scalar correction to the pressure
114 real energy = 0; //!< Correction to the energy
115 real dvdl = 0; //!< Correction to dH/dlambda
118 /*! \brief Computes and returns the dispersion correction for the pressure and energy
120 * \param[in] box The simulation unit cell
121 * \param[in] lambda The free-energy coupling parameter
123 Correction calculate(const matrix box, real lambda) const;
126 /*! \internal \brief Parameters that depend on the topology only
131 TopologyParams(const gmx_mtop_t& mtop,
132 const t_inputrec& inputrec,
135 gmx::ArrayRef<const real> nonbondedForceParameters);
137 //! The number of atoms for computing the atom density
138 int numAtomsForDensity_;
139 //! The number of interactions to correct for, usually num. atoms/2
140 real numCorrections_;
141 //! Average C6 coefficient for for topology A/B ([0]/[1])
142 std::array<real, 2> avcsix_;
143 //! Average C12 coefficient for for topology A/B ([0]/[1])
144 std::array<real, 2> avctwelve_;
147 /*! \internal \brief Parameters that depend on the interaction functions and topology
149 struct InteractionParams
152 ~InteractionParams();
154 //! Table used for correcting modified LJ interactions
155 std::unique_ptr<t_forcetable> dispersionCorrectionTable_;
157 //! Dispersion energy shift constant
158 real enershiftsix_ = 0;
159 //! Repulsion energy shift constant
160 real enershifttwelve_ = 0;
161 //! Dispersion energy difference per atom per unit of volume
162 real enerdiffsix_ = 0;
163 //! Repulsion energy difference per atom per unit of volume
164 real enerdifftwelve_ = 0;
165 //! Dispersion virial difference per atom per unit of volume
166 real virdiffsix_ = 0;
167 //! Repulsion virial difference per atom per unit of volume
168 real virdifftwelve_ = 0;
171 //! Sets the interaction parameters
172 static void setInteractionParameters(DispersionCorrection::InteractionParams* iParams,
173 const interaction_const_t& ic,
174 const char* tableFileName);
176 //! Returns whether we correct both dispersion and repulsion
177 bool correctFullInteraction() const;
179 //! Type of dispersion correction
181 //! Type of Van der Waals interaction
183 //! Free-energy perturbation
185 //! Topology parameters
186 TopologyParams topParams_;
187 //! Interaction parameters
188 InteractionParams iParams_;