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