Move foreign potential energy accumulation
[alexxy/gromacs.git] / src / gromacs / mdlib / enerdata_utils.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include "enerdata_utils.h"
41
42 #include "gromacs/mdtypes/enerdata.h"
43 #include "gromacs/mdtypes/inputrec.h"
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/utility/smalloc.h"
46
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)
52 {
53 }
54
55 static real sum_v(int n, gmx::ArrayRef<const real> v)
56 {
57     real t;
58     int  i;
59
60     t = 0.0;
61     for (i = 0; (i < n); i++)
62     {
63         t = t + v[i];
64     }
65
66     return t;
67 }
68
69 void sum_epot(const gmx_grppairener_t& grpp, real* epot)
70 {
71     int i;
72
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]);
78
79     /* lattice part of LR doesnt belong to any group
80      * and has been added earlier
81      */
82     epot[F_BHAM] = sum_v(grpp.nener, grpp.ener[egBHAMSR]);
83
84     epot[F_EPOT] = 0;
85     for (i = 0; (i < F_EPOT); i++)
86     {
87         if (i != F_DISRESVIOL && i != F_ORIRESDEV)
88         {
89             epot[F_EPOT] += epot[i];
90         }
91     }
92 }
93
94 /* Adds computed dV/dlambda contributions to the requested outputs
95  *
96  * Also adds the dispersion correction dV/dlambda to the VdW term.
97  * Note that kinetic energy and constraint contributions are handled in sum_dkdl()
98  */
99 static void sum_dvdl(gmx_enerdata_t* enerd, const t_lambda& fepvals)
100 {
101     // Add dispersion correction to the VdW component
102     enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW];
103
104     for (size_t i = 0; i < enerd->enerpart_lambda.size(); i++)
105     {
106         enerd->dhdlLambda[i] += enerd->term[F_DVDL_VDW];
107     }
108
109     enerd->term[F_DVDL] = 0.0;
110     for (int i = 0; i < efptNR; i++)
111     {
112         if (fepvals.separate_dvdl[i])
113         {
114             /* could this be done more readably/compactly? */
115             int index;
116             switch (i)
117             {
118                 case (efptMASS): index = F_DKDL; break;
119                 case (efptCOUL): index = F_DVDL_COUL; break;
120                 case (efptVDW): index = F_DVDL_VDW; break;
121                 case (efptBONDED): index = F_DVDL_BONDED; break;
122                 case (efptRESTRAINT): index = F_DVDL_RESTRAINT; break;
123                 default: index = F_DVDL; break;
124             }
125             enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
126             if (debug)
127             {
128                 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n", efpt_names[i], i,
129                         enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
130             }
131         }
132         else
133         {
134             enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
135             if (debug)
136             {
137                 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n", efpt_names[0], i,
138                         enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
139             }
140         }
141     }
142 }
143
144 void accumulatePotentialEnergies(gmx_enerdata_t* enerd, gmx::ArrayRef<const real> lambda, const t_lambda* fepvals)
145 {
146     sum_epot(enerd->grpp, enerd->term);
147
148     if (fepvals)
149     {
150         /* Note that sum_dvdl() adds the dispersion correction enerd->dvdl_lin[efptVDW],
151          * so sum_dvdl() should be called before computing the foreign lambda energy differences.
152          */
153         sum_dvdl(enerd, *fepvals);
154
155         /* Sum the foreign lambda energy difference contributions.
156          * Note that here we only add the potential energy components.
157          * The constraint and kinetic energy components are add after integration
158          * by sum_dhdl().
159          */
160         for (int i = 0; i < fepvals->n_lambda; i++)
161         {
162             /* note we are iterating over fepvals here!
163                For the current lam, dlam = 0 automatically,
164                so we don't need to add anything to the
165                enerd->enerpart_lambda[0] */
166
167             /* we don't need to worry about dvdl_lin contributions to dE at
168                current lambda, because the contributions to the current
169                lambda are automatically zeroed */
170
171             double& enerpart_lambda = enerd->enerpart_lambda[i + 1];
172
173             for (gmx::index j = 0; j < lambda.ssize(); j++)
174             {
175                 /* Note that this loop is over all dhdl components, not just the separated ones */
176                 const double dlam = fepvals->all_lambda[j][i] - lambda[j];
177
178                 enerpart_lambda += dlam * enerd->dvdl_lin[j];
179             }
180         }
181     }
182 }
183
184 void accumulateKineticLambdaComponents(gmx_enerdata_t*           enerd,
185                                        gmx::ArrayRef<const real> lambda,
186                                        const t_lambda&           fepvals)
187 {
188     if (fepvals.separate_dvdl[efptBONDED])
189     {
190         enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
191     }
192     else
193     {
194         enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
195     }
196
197     for (int i = 0; i < fepvals.n_lambda; i++)
198     {
199         /* note we are iterating over fepvals here!
200            For the current lam, dlam = 0 automatically,
201            so we don't need to add anything to the
202            enerd->enerpart_lambda[0] */
203
204         double& enerpart_lambda = enerd->enerpart_lambda[i + 1];
205
206         /* Note that potential energy terms have been added by sum_epot() -> sum_dvdl() */
207
208         /* Constraints can not be evaluated at foreign lambdas, so we add
209          * a linear extrapolation. This is an approximation, but usually
210          * quite accurate since constraints change little between lambdas.
211          */
212         const int    lambdaIndex = (fepvals.separate_dvdl[efptBONDED] ? efptBONDED : efptFEP);
213         const double dlam        = fepvals.all_lambda[lambdaIndex][i] - lambda[lambdaIndex];
214         enerpart_lambda += dlam * enerd->term[F_DVDL_CONSTR];
215
216         if (!fepvals.separate_dvdl[efptMASS])
217         {
218             const double dlam = fepvals.all_lambda[efptMASS][i] - lambda[efptMASS];
219             enerpart_lambda += dlam * enerd->term[F_DKDL];
220         }
221     }
222
223     /* The constrain contribution is now included in other terms, so clear it */
224     enerd->term[F_DVDL_CONSTR] = 0;
225 }
226
227
228 void reset_foreign_enerdata(gmx_enerdata_t* enerd)
229 {
230     int i, j;
231
232     /* First reset all foreign energy components.  Foreign energies always called on
233        neighbor search steps */
234     for (i = 0; (i < egNR); i++)
235     {
236         for (j = 0; (j < enerd->grpp.nener); j++)
237         {
238             enerd->foreign_grpp.ener[i][j] = 0.0;
239         }
240     }
241
242     /* potential energy components */
243     for (i = 0; (i <= F_EPOT); i++)
244     {
245         enerd->foreign_term[i] = 0.0;
246     }
247 }
248
249 void reset_dvdl_enerdata(gmx_enerdata_t* enerd)
250 {
251     for (int i = 0; i < efptNR; i++)
252     {
253         enerd->dvdl_lin[i]    = 0.0;
254         enerd->dvdl_nonlin[i] = 0.0;
255     }
256 }
257
258 void reset_enerdata(gmx_enerdata_t* enerd)
259 {
260     int i, j;
261
262     /* First reset all energy components. */
263     for (i = 0; (i < egNR); i++)
264     {
265         for (j = 0; (j < enerd->grpp.nener); j++)
266         {
267             enerd->grpp.ener[i][j] = 0.0_real;
268         }
269     }
270
271     /* Normal potential energy components */
272     for (i = 0; (i <= F_EPOT); i++)
273     {
274         enerd->term[i] = 0.0_real;
275     }
276     enerd->term[F_PDISPCORR]      = 0.0_real;
277     enerd->term[F_DVDL]           = 0.0_real;
278     enerd->term[F_DVDL_COUL]      = 0.0_real;
279     enerd->term[F_DVDL_VDW]       = 0.0_real;
280     enerd->term[F_DVDL_BONDED]    = 0.0_real;
281     enerd->term[F_DVDL_RESTRAINT] = 0.0_real;
282     enerd->term[F_DKDL]           = 0.0_real;
283     std::fill(enerd->enerpart_lambda.begin(), enerd->enerpart_lambda.end(), 0);
284     std::fill(enerd->dhdlLambda.begin(), enerd->dhdlLambda.end(), 0);
285     /* reset foreign energy data and dvdl - separate functions since they are also called elsewhere */
286     reset_foreign_enerdata(enerd);
287     reset_dvdl_enerdata(enerd);
288 }