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