2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020,2021, 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;
50 enum class DispersionCorrectionType : int;
51 enum class VanDerWaalsType : int;
52 enum class FreeEnergyPerturbationType : int;
60 class DispersionCorrection
63 /*! \brief Constructor
65 * \param[in] mtop The global topology
66 * \param[in] inputrec The input record
67 * \param[in] useBuckingham True when Buckingham is used instead of LJ
68 * \param[in] numAtomTypes The number of non-bonded atom types
69 * \param[in] nonbondedForceParameters The LJ or Bham parameter matrix stored as a flat list
70 * \param[in] ic The nonbonded interaction parameters
71 * \param[in] tableFileName Table file name, should != nullptr (checked)
73 DispersionCorrection(const gmx_mtop_t& mtop,
74 const t_inputrec& inputrec,
77 gmx::ArrayRef<const real> nonbondedForceParameters,
78 const interaction_const_t& ic,
79 const char* tableFileName);
81 /*! \brief Print dispersion correction information to the log file
83 * \param[in] mdlog The MD logger
85 void print(const gmx::MDLogger& mdlog) const;
87 /*! \brief Computes and sets energy and virial correction parameters
89 * Sets all parameters that are affected by the cut-off and/or the
90 * LJ-Ewald coefficient. Should be called before calling calculate()
91 * and whenever interaction settings change, e.g. PME tuning.
93 * \param[in] ic The nonbonded interaction parameters
95 void setParameters(const interaction_const_t& ic);
98 * \brief Struct for returning all dispersion correction quantities
102 /*! \brief Correct the virial tensor for the missing dispersion
104 * \param[in,out] virialTensor The virial tensor to correct
106 void correctVirial(tensor virialTensor) const
108 for (int m = 0; m < DIM; m++)
110 virialTensor[m][m] += virial;
114 real virial = 0; //!< Scalar correction to the virial
115 real pressure = 0; //!< Scalar correction to the pressure
116 real energy = 0; //!< Correction to the energy
117 real dvdl = 0; //!< Correction to dH/dlambda
120 /*! \brief Computes and returns the dispersion correction for the pressure and energy
122 * \param[in] box The simulation unit cell
123 * \param[in] lambda The free-energy coupling parameter
125 Correction calculate(const matrix box, real lambda) const;
128 /*! \internal \brief Parameters that depend on the topology only
133 TopologyParams(const gmx_mtop_t& mtop,
134 const t_inputrec& inputrec,
137 gmx::ArrayRef<const real> nonbondedForceParameters);
139 //! The number of atoms for computing the atom density
140 int numAtomsForDensity_;
141 //! The number of interactions to correct for, usually num. atoms/2
142 real numCorrections_;
143 //! Average C6 coefficient for for topology A/B ([0]/[1])
144 std::array<real, 2> avcsix_;
145 //! Average C12 coefficient for for topology A/B ([0]/[1])
146 std::array<real, 2> avctwelve_;
149 /*! \internal \brief Parameters that depend on the interaction functions and topology
151 struct InteractionParams
154 ~InteractionParams();
156 //! Table used for correcting modified LJ interactions
157 std::unique_ptr<t_forcetable> dispersionCorrectionTable_;
159 //! Dispersion energy shift constant
160 real enershiftsix_ = 0;
161 //! Repulsion energy shift constant
162 real enershifttwelve_ = 0;
163 //! Dispersion energy difference per atom per unit of volume
164 real enerdiffsix_ = 0;
165 //! Repulsion energy difference per atom per unit of volume
166 real enerdifftwelve_ = 0;
167 //! Dispersion virial difference per atom per unit of volume
168 real virdiffsix_ = 0;
169 //! Repulsion virial difference per atom per unit of volume
170 real virdifftwelve_ = 0;
173 //! Sets the interaction parameters
174 static void setInteractionParameters(DispersionCorrection::InteractionParams* iParams,
175 const interaction_const_t& ic,
176 const char* tableFileName);
178 //! Returns whether we correct both dispersion and repulsion
179 bool correctFullInteraction() const;
181 //! Type of dispersion correction
182 DispersionCorrectionType eDispCorr_;
183 //! Type of Van der Waals interaction
184 VanDerWaalsType vdwType_;
185 //! Free-energy perturbation
186 FreeEnergyPerturbationType eFep_;
187 //! Topology parameters
188 TopologyParams topParams_;
189 //! Interaction parameters
190 InteractionParams iParams_;