88c8c7b40b9c2ed409a13711e56d2d29ed8d94b1
[alexxy/gromacs.git] / src / gromacs / mdlib / dispersioncorrection.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #ifndef GMX_MDLIB_DISPERSIONCORRECTION_H
36 #define GMX_MDLIB_DISPERSIONCORRECTION_H
37
38 #include <cstdio>
39
40 #include <memory>
41
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/utility/arrayref.h"
44
45 struct gmx_mtop_t;
46 struct interaction_const_t;
47 struct t_forcerec;
48 struct t_forcetable;
49 struct t_inputrec;
50
51 namespace gmx
52 {
53 class MDLogger;
54 } // namespace gmx
55
56 class DispersionCorrection
57 {
58 public:
59     /*! \brief Constructor
60      *
61      * \param[in] mtop           The global topology
62      * \param[in] inputrec       The input record
63      * \param[in] useBuckingham  True when Buckingham is used instead of LJ
64      * \param[in] numAtomTypes   The number of non-bonded atom types
65      * \param[in] nonbondedForceParameters  The LJ or Bham parameter matrix stored as a flat list
66      * \param[in] ic             The nonbonded interaction parameters
67      * \param[in] tableFileName  Table file name, should != nullptr (checked)
68      */
69     DispersionCorrection(const gmx_mtop_t&          mtop,
70                          const t_inputrec&          inputrec,
71                          bool                       useBuckingham,
72                          int                        numAtomTypes,
73                          gmx::ArrayRef<const real>  nonbondedForceParameters,
74                          const interaction_const_t& ic,
75                          const char*                tableFileName);
76
77     /*! \brief Print dispersion correction information to the log file
78      *
79      * \param[in] mdlog  The MD logger
80      */
81     void print(const gmx::MDLogger& mdlog) const;
82
83     /*! \brief Computes and sets energy and virial correction parameters
84      *
85      * Sets all parameters that are affected by the cut-off and/or the
86      * LJ-Ewald coefficient. Should be called before calling calculate()
87      * and whenever interaction settings change, e.g. PME tuning.
88      *
89      * \param[in] ic  The nonbonded interaction parameters
90      */
91     void setParameters(const interaction_const_t& ic);
92
93     /*! \internal
94      * \brief Struct for returning all dispersion correction quantities
95      */
96     struct Correction
97     {
98         /*! \brief Correct the virial tensor for the missing dispersion
99          *
100          * \param[in,out] virialTensor  The virial tensor to correct
101          */
102         void correctVirial(tensor virialTensor) const
103         {
104             for (int m = 0; m < DIM; m++)
105             {
106                 virialTensor[m][m] += virial;
107             }
108         }
109
110         /*! \brief Correct the pressure tensor for the missing dispersion
111          *
112          * \param[in,out] pressureTensor  The pressure tensor to correct
113          */
114         void correctPressure(tensor pressureTensor) const
115         {
116             for (int m = 0; m < DIM; m++)
117             {
118                 pressureTensor[m][m] += pressure;
119             }
120         }
121
122         real virial   = 0; //!< Scalar correction to the virial
123         real pressure = 0; //!< Scalar correction to the pressure
124         real energy   = 0; //!< Correction to the energy
125         real dvdl     = 0; //!< Correction to dH/dlambda
126     };
127
128     /*! \brief Computes and returns the dispersion correction for the pressure and energy
129      *
130      * \param[in]  box       The simulation unit cell
131      * \param[in]  lambda    The free-energy coupling parameter
132      */
133     Correction calculate(const matrix box, real lambda) const;
134
135 private:
136     /*! \internal \brief Parameters that depend on the topology only
137      */
138     class TopologyParams
139     {
140     public:
141         TopologyParams(const gmx_mtop_t&         mtop,
142                        const t_inputrec&         inputrec,
143                        bool                      useBuckingham,
144                        int                       numAtomTypes,
145                        gmx::ArrayRef<const real> nonbondedForceParameters);
146
147         //! The number of atoms for computing the atom density
148         int numAtomsForDensity_;
149         //! The number of interactions to correct for, usually num. atoms/2
150         real numCorrections_;
151         //! Average C6 coefficient for for topology A/B ([0]/[1])
152         std::array<real, 2> avcsix_;
153         //! Average C12 coefficient for for topology A/B ([0]/[1])
154         std::array<real, 2> avctwelve_;
155     };
156
157     /*! \internal \brief Parameters that depend on the interaction functions and topology
158      */
159     struct InteractionParams
160     {
161     public:
162         ~InteractionParams();
163
164         //! Table used for correcting modified LJ interactions
165         std::unique_ptr<t_forcetable> dispersionCorrectionTable_;
166
167         //! Dispersion energy shift constant
168         real enershiftsix_ = 0;
169         //! Repulsion energy shift constant
170         real enershifttwelve_ = 0;
171         //! Dispersion energy difference per atom per unit of volume
172         real enerdiffsix_ = 0;
173         //! Repulsion energy difference per atom per unit of volume
174         real enerdifftwelve_ = 0;
175         //! Dispersion virial difference per atom per unit of volume
176         real virdiffsix_ = 0;
177         //! Repulsion virial difference per atom per unit of volume
178         real virdifftwelve_ = 0;
179     };
180
181     //! Sets the interaction parameters
182     static void setInteractionParameters(DispersionCorrection::InteractionParams* iParams,
183                                          const interaction_const_t&               ic,
184                                          const char*                              tableFileName);
185
186     //! Returns whether we correct both dispersion and repulsion
187     bool correctFullInteraction() const;
188
189     //! Type of dispersion correction
190     int eDispCorr_;
191     //! Type of Van der Waals interaction
192     int vdwType_;
193     //! Free-energy perturbation
194     int eFep_;
195     //! Topology parameters
196     TopologyParams topParams_;
197     //! Interaction parameters
198     InteractionParams iParams_;
199 };
200
201 #endif