Make non bonded energy terms enum class
[alexxy/gromacs.git] / src / gromacs / mdtypes / enerdata.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,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.
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_MDTYPES_TYPES_ENERDATA_H
36 #define GMX_MDTYPES_TYPES_ENERDATA_H
37
38 #include <array>
39 #include <utility>
40 #include <vector>
41
42 #include "gromacs/mdtypes/md_enums.h"
43 #include "gromacs/topology/idef.h"
44 #include "gromacs/utility/arrayref.h"
45 #include "gromacs/utility/enumerationhelpers.h"
46 #include "gromacs/utility/gmxassert.h"
47 #include "gromacs/utility/real.h"
48
49 struct t_commrec;
50 struct t_lambda;
51
52 // The non-bonded energy terms accumulated for energy group pairs
53 enum class NonBondedEnergyTerms : int
54 {
55     CoulombSR,
56     LJSR,
57     BuckinghamSR,
58     Coulomb14,
59     LJ14,
60     Count
61 };
62
63 // Struct for accumulating non-bonded energies between energy group pairs
64 struct gmx_grppairener_t
65 {
66     gmx_grppairener_t(int numEnergyGroups) : nener(numEnergyGroups * numEnergyGroups)
67     {
68         for (auto& term : energyGroupPairTerms)
69         {
70             term.resize(nener);
71         }
72     }
73
74     int nener; /* The number of energy group pairs */
75     gmx::EnumerationArray<NonBondedEnergyTerms, std::vector<real>> energyGroupPairTerms; /* Energy terms for each pair of groups */
76 };
77
78 //! Accumulates free-energy foreign lambda energies and dH/dlamba
79 class ForeignLambdaTerms
80 {
81 public:
82     /*! \brief Constructor
83      *
84      * \param[in] numLambdas  The number of foreign lambda values
85      */
86     ForeignLambdaTerms(int numLambdas);
87
88     //! Returns the number of foreign lambda values
89     int numLambdas() const { return numLambdas_; }
90
91     //! Returns the H(lambdaIndex) - H(lambda_current)
92     double deltaH(int lambdaIndex) const { return energies_[1 + lambdaIndex] - energies_[0]; }
93
94     /*! \brief Returns a list of partial energies, the part which depends on lambda),
95      * current lambda in entry 0, foreign lambda i in entry 1+i
96      *
97      * Note: the potential terms needs to be finalized before calling this method.
98      */
99     gmx::ArrayRef<double> energies()
100     {
101         GMX_ASSERT(finalizedPotentialContributions_, "Should be finalized");
102         return energies_;
103     }
104
105     /*! \brief Returns a list of partial energies, the part which depends on lambda),
106      * current lambda in entry 0, foreign lambda i in entry 1+i
107      *
108      * Note: the potential terms needs to be finalized before calling this method.
109      */
110     gmx::ArrayRef<const double> energies() const
111     {
112         GMX_ASSERT(finalizedPotentialContributions_, "Should be finalized");
113         return energies_;
114     }
115
116     /*! \brief Adds an energy and dV/dl constribution to lambda list index \p listIndex
117      *
118      * This should only be used for terms with non-linear dependence on lambda
119      * The value passed as listIndex should be 0 for the current lambda
120      * and 1+i for foreign lambda index i.
121      */
122     void accumulate(int listIndex, double energy, double dvdl)
123     {
124         GMX_ASSERT(!finalizedPotentialContributions_,
125                    "Can only accumulate with an unfinalized object");
126
127         energies_[listIndex] += energy;
128         dhdl_[listIndex] += dvdl;
129     }
130
131     /*! \brief Finalizes the potential (non-kinetic) terms
132      *
133      * Note: This can be called multiple times during the same force calculations
134      * without affecting the results.
135      *
136      * \param[in] dvdlLinear  List of dV/dlambda contributions of size efptNR with depend linearly on lambda
137      * \param[in] lambda      Lambda values for the efptNR contribution types
138      * \param[in] fepvals     Free-energy parameters
139      */
140     void finalizePotentialContributions(gmx::ArrayRef<const double> dvdlLinear,
141                                         gmx::ArrayRef<const real>   lambda,
142                                         const t_lambda&             fepvals);
143
144     /*! \brief Accumulates the kinetic and constraint free-energy contributions
145      *
146      * \param[in] energyTerms  List of energy terms, pass \p term in \p gmx_enerdata_t
147      * \param[in] dhdlMass     The mass dependent contribution to dH/dlambda
148      * \param[in] lambda       Lambda values for the efptNR contribution types
149      * \param[in] fepvals      Free-energy parameters
150      */
151     void finalizeKineticContributions(gmx::ArrayRef<const real> energyTerms,
152                                       double                    dhdlMass,
153                                       gmx::ArrayRef<const real> lambda,
154                                       const t_lambda&           fepvals);
155
156     /*! \brief Returns a pair of lists of deltaH and dH/dlambda
157      *
158      * Both lists are of size numLambdas() and are indexed with the lambda index.
159      *
160      * Note: should only be called after the object has been finalized by a call to
161      * accumulateLinearPotentialComponents() (is asserted).
162      *
163      * \param[in] cr  Communication record, used to reduce the terms when !=nullptr
164      */
165     std::pair<std::vector<double>, std::vector<double>> getTerms(const t_commrec* cr) const;
166
167     //! Sets all terms to 0
168     void zeroAllTerms();
169
170 private:
171     //! As accumulate(), but for kinetic contributions
172     void accumulateKinetic(int listIndex, double energy, double dhdl);
173
174     //! Add a dH/dl contribution that does not depend on lambda to all foreign dH/dl terms
175     void addConstantDhdl(double dhdl);
176
177     //! The number of foreign lambdas
178     int numLambdas_;
179     //! Storage for foreign lambda energies
180     std::vector<double> energies_;
181     //! Storage for foreign lambda dH/dlambda
182     std::vector<double> dhdl_;
183     //! Tells whether all potential energy contributions have been accumulated
184     bool finalizedPotentialContributions_ = false;
185 };
186
187 //! Struct for accumulating all potential energy terms and some kinetic energy terms
188 struct gmx_enerdata_t
189 {
190     /*! \brief
191      * Constructor with specific number of energy groups and lambdas.
192      *
193      * \param[in] numEnergyGroups Number of energy groups used.
194      * \param[in] numFepLambdas   Number of free energy lambdas, zero if none.
195      */
196     gmx_enerdata_t(int numEnergyGroups, int numFepLambdas);
197
198     //! The energies for all different interaction types
199     real term[F_NRE] = { 0 };
200     //! Energy group pair non-bonded energies
201     struct gmx_grppairener_t grpp;
202     //! Contributions to dV/dlambda with linear dependence on lambda
203     gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, double> dvdl_lin = { 0 };
204     //! Contributions to dV/dlambda with non-linear dependence on lambda
205     gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, double> dvdl_nonlin = { 0 };
206     /* The idea is that dvdl terms with linear lambda dependence will be added
207      * automatically to enerpart_lambda. Terms with non-linear lambda dependence
208      * should explicitly determine the energies at foreign lambda points
209      * when n_lambda > 0. */
210
211     //! Foreign lambda energies and dH/dl
212     ForeignLambdaTerms foreignLambdaTerms;
213
214     //! Alternate, temporary array for storing foreign lambda energies
215     real foreign_term[F_NRE] = { 0 };
216     //! Alternate, temporary  array for storing foreign lambda group pair energies
217     struct gmx_grppairener_t foreign_grpp;
218 };
219
220 #endif