2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 #include "enerdata_utils.h"
42 #include "gromacs/gmxlib/network.h"
43 #include "gromacs/mdtypes/commrec.h"
44 #include "gromacs/mdtypes/enerdata.h"
45 #include "gromacs/mdtypes/inputrec.h"
46 #include "gromacs/utility/enumerationhelpers.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/smalloc.h"
50 ForeignLambdaTerms::ForeignLambdaTerms(int numLambdas) :
51 numLambdas_(numLambdas), energies_(1 + numLambdas), dhdl_(1 + numLambdas)
55 std::pair<std::vector<double>, std::vector<double>> ForeignLambdaTerms::getTerms(const t_commrec* cr) const
57 GMX_RELEASE_ASSERT(finalizedPotentialContributions_,
58 "The object needs to be finalized before calling getTerms");
60 std::vector<double> data(2 * numLambdas_);
61 for (int i = 0; i < numLambdas_; i++)
63 data[i] = energies_[1 + i] - energies_[0];
64 data[numLambdas_ + i] = dhdl_[1 + i];
66 if (cr && cr->nnodes > 1)
68 gmx_sumd(data.size(), data.data(), cr);
70 auto dataMid = data.begin() + numLambdas_;
72 return { { data.begin(), dataMid }, { dataMid, data.end() } };
75 void ForeignLambdaTerms::zeroAllTerms()
77 std::fill(energies_.begin(), energies_.end(), 0.0);
78 std::fill(dhdl_.begin(), dhdl_.end(), 0.0);
79 finalizedPotentialContributions_ = false;
82 gmx_enerdata_t::gmx_enerdata_t(int numEnergyGroups, int numFepLambdas) :
83 grpp(numEnergyGroups), foreignLambdaTerms(numFepLambdas)
87 static real sum_v(int n, gmx::ArrayRef<const real> v)
93 for (i = 0; (i < n); i++)
101 void sum_epot(const gmx_grppairener_t& grpp, real* epot)
105 /* Accumulate energies */
106 epot[F_COUL_SR] = sum_v(grpp.nener, grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR]);
107 epot[F_LJ] = sum_v(grpp.nener, grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR]);
108 epot[F_LJ14] = sum_v(grpp.nener, grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJ14]);
109 epot[F_COUL14] = sum_v(grpp.nener, grpp.energyGroupPairTerms[NonBondedEnergyTerms::Coulomb14]);
111 /* lattice part of LR doesnt belong to any group
112 * and has been added earlier
114 epot[F_BHAM] = sum_v(grpp.nener, grpp.energyGroupPairTerms[NonBondedEnergyTerms::BuckinghamSR]);
117 for (i = 0; (i < F_EPOT); i++)
119 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
121 epot[F_EPOT] += epot[i];
126 // Adds computed dH/dlambda contribution i to the requested output
127 static void set_dhdl_output(gmx_enerdata_t* enerd, FreeEnergyPerturbationCouplingType i, const t_lambda& fepvals)
129 if (fepvals.separate_dvdl[i])
131 /* Translate free-energy term indices to idef term indices.
132 * Could this be done more readably/compactly?
137 case (FreeEnergyPerturbationCouplingType::Mass): index = F_DKDL; break;
138 case (FreeEnergyPerturbationCouplingType::Coul): index = F_DVDL_COUL; break;
139 case (FreeEnergyPerturbationCouplingType::Vdw): index = F_DVDL_VDW; break;
140 case (FreeEnergyPerturbationCouplingType::Bonded): index = F_DVDL_BONDED; break;
141 case (FreeEnergyPerturbationCouplingType::Restraint): index = F_DVDL_RESTRAINT; break;
142 default: index = F_DVDL; break;
144 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
148 "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n",
149 enumValueToString(i),
152 enerd->dvdl_nonlin[i],
158 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
162 "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n",
163 enumValueToString(FreeEnergyPerturbationCouplingType::Fep),
166 enerd->dvdl_nonlin[i],
172 void ForeignLambdaTerms::addConstantDhdl(const double dhdl)
174 for (double& foreignDhdl : dhdl_)
180 void ForeignLambdaTerms::finalizePotentialContributions(gmx::ArrayRef<const double> dvdlLinear,
181 gmx::ArrayRef<const real> lambda,
182 const t_lambda& fepvals)
184 if (finalizedPotentialContributions_)
190 for (int i = 0; i < static_cast<int>(FreeEnergyPerturbationCouplingType::Count); i++)
192 dvdl_lin += dvdlLinear[i];
194 addConstantDhdl(dvdl_lin);
196 /* Sum the foreign lambda energy difference contributions.
197 * Note that here we only add the potential energy components.
198 * The constraint and kinetic energy components are add after integration
201 for (int i = 0; i < fepvals.n_lambda; i++)
203 /* note we are iterating over fepvals here!
204 For the current lam, dlam = 0 automatically,
205 so we don't need to add anything to the
206 enerd->enerpart_lambda[0] */
208 /* we don't need to worry about dvdl_lin contributions to dE at
209 current lambda, because the contributions to the current
210 lambda are automatically zeroed */
212 double enerpart_lambda = 0;
213 for (gmx::index j = 0; j < lambda.ssize(); j++)
215 /* Note that this loop is over all dhdl components, not just the separated ones */
216 const double dlam = fepvals.all_lambda[j][i] - lambda[j];
218 enerpart_lambda += dlam * dvdlLinear[j];
220 accumulate(1 + i, enerpart_lambda, 0);
223 finalizedPotentialContributions_ = true;
226 void accumulatePotentialEnergies(gmx_enerdata_t* enerd, gmx::ArrayRef<const real> lambda, const t_lambda* fepvals)
228 sum_epot(enerd->grpp, enerd->term.data());
232 enerd->term[F_DVDL] = 0.0;
233 for (auto i : gmx::EnumerationWrapper<FreeEnergyPerturbationCouplingType>{})
235 // Skip kinetic terms here, as those are not available here yet
236 if (i != FreeEnergyPerturbationCouplingType::Mass)
238 set_dhdl_output(enerd, i, *fepvals);
242 enerd->foreignLambdaTerms.finalizePotentialContributions(enerd->dvdl_lin, lambda, *fepvals);
246 void ForeignLambdaTerms::accumulateKinetic(int listIndex, double energy, double dhdl)
248 energies_[listIndex] += energy;
249 dhdl_[listIndex] += dhdl;
252 void ForeignLambdaTerms::finalizeKineticContributions(gmx::ArrayRef<const real> energyTerms,
253 const double dhdlMass,
254 gmx::ArrayRef<const real> lambda,
255 const t_lambda& fepvals)
257 // Add perturbed mass contributions
258 addConstantDhdl(dhdlMass);
260 // Treat current lambda, the deltaH contribution is 0 as delta-lambda=0 for the current lambda
261 accumulateKinetic(0, 0.0, energyTerms[F_DVDL_CONSTR]);
262 if (!fepvals.separate_dvdl[FreeEnergyPerturbationCouplingType::Mass])
264 accumulateKinetic(0, 0.0, energyTerms[F_DKDL]);
267 for (int i = 0; i < fepvals.n_lambda; i++)
269 /* Note that potential energy terms have been added by sum_epot() -> sum_dvdl() */
271 /* Constraints can not be evaluated at foreign lambdas, so we add
272 * a linear extrapolation. This is an approximation, but usually
273 * quite accurate since constraints change little between lambdas.
275 const FreeEnergyPerturbationCouplingType lambdaIndex =
276 (fepvals.separate_dvdl[FreeEnergyPerturbationCouplingType::Bonded]
277 ? FreeEnergyPerturbationCouplingType::Bonded
278 : FreeEnergyPerturbationCouplingType::Fep);
279 const double dlam = fepvals.all_lambda[lambdaIndex][i] - lambda[static_cast<int>(lambdaIndex)];
280 accumulateKinetic(1 + i, dlam * energyTerms[F_DVDL_CONSTR], energyTerms[F_DVDL_CONSTR]);
282 if (!fepvals.separate_dvdl[FreeEnergyPerturbationCouplingType::Mass])
284 const double dlam = fepvals.all_lambda[FreeEnergyPerturbationCouplingType::Mass][i]
285 - lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Mass)];
286 accumulateKinetic(1 + i, dlam * energyTerms[F_DKDL], energyTerms[F_DKDL]);
291 void accumulateKineticLambdaComponents(gmx_enerdata_t* enerd,
292 gmx::ArrayRef<const real> lambda,
293 const t_lambda& fepvals)
295 if (fepvals.separate_dvdl[FreeEnergyPerturbationCouplingType::Bonded])
297 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
301 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
304 // Add computed mass dH/dlambda contribution to the requested output
305 set_dhdl_output(enerd, FreeEnergyPerturbationCouplingType::Mass, fepvals);
307 enerd->foreignLambdaTerms.finalizeKineticContributions(
308 enerd->term, enerd->dvdl_lin[FreeEnergyPerturbationCouplingType::Mass], lambda, fepvals);
310 /* The constrain contribution is now included in other terms, so clear it */
311 enerd->term[F_DVDL_CONSTR] = 0;
314 void gmx_grppairener_t::clear()
316 for (int i = 0; (i < static_cast<int>(NonBondedEnergyTerms::Count)); i++)
318 for (int j = 0; (j < nener); j++)
320 energyGroupPairTerms[i][j] = 0.0;
325 void reset_dvdl_enerdata(gmx_enerdata_t* enerd)
327 for (auto i : keysOf(enerd->dvdl_lin))
329 enerd->dvdl_lin[i] = 0.0;
330 enerd->dvdl_nonlin[i] = 0.0;
334 void reset_enerdata(gmx_enerdata_t* enerd)
338 /* First reset all energy components. */
339 for (i = 0; (i < static_cast<int>(NonBondedEnergyTerms::Count)); i++)
341 for (j = 0; (j < enerd->grpp.nener); j++)
343 enerd->grpp.energyGroupPairTerms[i][j] = 0.0_real;
347 /* Normal potential energy components */
348 for (i = 0; (i <= F_EPOT); i++)
350 enerd->term[i] = 0.0_real;
352 enerd->term[F_PDISPCORR] = 0.0_real;
353 enerd->term[F_DVDL] = 0.0_real;
354 enerd->term[F_DVDL_COUL] = 0.0_real;
355 enerd->term[F_DVDL_VDW] = 0.0_real;
356 enerd->term[F_DVDL_BONDED] = 0.0_real;
357 enerd->term[F_DVDL_RESTRAINT] = 0.0_real;
358 enerd->term[F_DKDL] = 0.0_real;
359 enerd->foreignLambdaTerms.zeroAllTerms();
360 /* reset dvdl - separate function since it is also called elsewhere */
361 reset_dvdl_enerdata(enerd);