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, 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/mdtypes/enerdata.h"
43 #include "gromacs/mdtypes/inputrec.h"
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/utility/smalloc.h"
47 gmx_enerdata_t::gmx_enerdata_t(int numEnergyGroups, int numFepLambdas) :
48 grpp(numEnergyGroups),
49 enerpart_lambda(numFepLambdas == 0 ? 0 : numFepLambdas + 1),
50 dhdlLambda(numFepLambdas == 0 ? 0 : numFepLambdas + 1),
51 foreign_grpp(numEnergyGroups)
55 static real sum_v(int n, gmx::ArrayRef<const real> v)
61 for (i = 0; (i < n); i++)
69 void sum_epot(gmx_grppairener_t* grpp, real* epot)
73 /* Accumulate energies */
74 epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
75 epot[F_LJ] = sum_v(grpp->nener, grpp->ener[egLJSR]);
76 epot[F_LJ14] = sum_v(grpp->nener, grpp->ener[egLJ14]);
77 epot[F_COUL14] = sum_v(grpp->nener, grpp->ener[egCOUL14]);
79 /* lattice part of LR doesnt belong to any group
80 * and has been added earlier
82 epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
85 for (i = 0; (i < F_EPOT); i++)
87 if (i != F_DISRESVIOL && i != F_ORIRESDEV)
89 epot[F_EPOT] += epot[i];
94 void sum_dhdl(gmx_enerdata_t* enerd, gmx::ArrayRef<const real> lambda, const t_lambda& fepvals)
98 enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
100 for (size_t i = 0; i < enerd->enerpart_lambda.size(); i++)
102 enerd->dhdlLambda[i] += enerd->term[F_DVDL_VDW];
105 enerd->term[F_DVDL] = 0.0;
106 for (int i = 0; i < efptNR; i++)
108 if (fepvals.separate_dvdl[i])
110 /* could this be done more readably/compactly? */
113 case (efptMASS): index = F_DKDL; break;
114 case (efptCOUL): index = F_DVDL_COUL; break;
115 case (efptVDW): index = F_DVDL_VDW; break;
116 case (efptBONDED): index = F_DVDL_BONDED; break;
117 case (efptRESTRAINT): index = F_DVDL_RESTRAINT; break;
118 default: index = F_DVDL; break;
120 enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
123 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n", efpt_names[i], i,
124 enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
129 enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
132 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n", efpt_names[0], i,
133 enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
138 if (fepvals.separate_dvdl[efptBONDED])
140 enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
144 enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
147 for (int i = 0; i < fepvals.n_lambda; i++)
149 /* note we are iterating over fepvals here!
150 For the current lam, dlam = 0 automatically,
151 so we don't need to add anything to the
152 enerd->enerpart_lambda[0] */
154 /* we don't need to worry about dvdl_lin contributions to dE at
155 current lambda, because the contributions to the current
156 lambda are automatically zeroed */
158 double& enerpart_lambda = enerd->enerpart_lambda[i + 1];
160 for (gmx::index j = 0; j < lambda.ssize(); j++)
162 /* Note that this loop is over all dhdl components, not just the separated ones */
163 const double dlam = fepvals.all_lambda[j][i] - lambda[j];
165 enerpart_lambda += dlam * enerd->dvdl_lin[j];
167 /* Constraints can not be evaluated at foreign lambdas, so we add
168 * a linear extrapolation. This is an approximation, but usually
169 * quite accurate since constraints change little between lambdas.
171 if ((j == efptBONDED && fepvals.separate_dvdl[efptBONDED])
172 || (j == efptFEP && !fepvals.separate_dvdl[efptBONDED]))
174 enerpart_lambda += dlam * enerd->term[F_DVDL_CONSTR];
177 if (j == efptMASS && !fepvals.separate_dvdl[j])
179 enerpart_lambda += dlam * enerd->term[F_DKDL];
184 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
185 fepvals.all_lambda[j][i], efpt_names[j],
186 enerpart_lambda - enerd->enerpart_lambda[0], dlam, enerd->dvdl_lin[j]);
191 /* The constrain contribution is now included in other terms, so clear it */
192 enerd->term[F_DVDL_CONSTR] = 0;
196 void reset_foreign_enerdata(gmx_enerdata_t* enerd)
200 /* First reset all foreign energy components. Foreign energies always called on
201 neighbor search steps */
202 for (i = 0; (i < egNR); i++)
204 for (j = 0; (j < enerd->grpp.nener); j++)
206 enerd->foreign_grpp.ener[i][j] = 0.0;
210 /* potential energy components */
211 for (i = 0; (i <= F_EPOT); i++)
213 enerd->foreign_term[i] = 0.0;
217 void reset_dvdl_enerdata(gmx_enerdata_t* enerd)
219 for (int i = 0; i < efptNR; i++)
221 enerd->dvdl_lin[i] = 0.0;
222 enerd->dvdl_nonlin[i] = 0.0;
226 void reset_enerdata(gmx_enerdata_t* enerd)
230 /* First reset all energy components. */
231 for (i = 0; (i < egNR); i++)
233 for (j = 0; (j < enerd->grpp.nener); j++)
235 enerd->grpp.ener[i][j] = 0.0_real;
239 /* Normal potential energy components */
240 for (i = 0; (i <= F_EPOT); i++)
242 enerd->term[i] = 0.0_real;
244 enerd->term[F_PDISPCORR] = 0.0_real;
245 enerd->term[F_DVDL] = 0.0_real;
246 enerd->term[F_DVDL_COUL] = 0.0_real;
247 enerd->term[F_DVDL_VDW] = 0.0_real;
248 enerd->term[F_DVDL_BONDED] = 0.0_real;
249 enerd->term[F_DVDL_RESTRAINT] = 0.0_real;
250 enerd->term[F_DKDL] = 0.0_real;
251 std::fill(enerd->enerpart_lambda.begin(), enerd->enerpart_lambda.end(), 0);
252 std::fill(enerd->dhdlLambda.begin(), enerd->dhdlLambda.end(), 0);
253 /* reset foreign energy data and dvdl - separate functions since they are also called elsewhere */
254 reset_foreign_enerdata(enerd);
255 reset_dvdl_enerdata(enerd);