ee70f7a1301d88bfb2a0995cd51dd3e4cb93dd9f
[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,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "enerdata_utils.h"
40
41 #include "gromacs/mdtypes/enerdata.h"
42 #include "gromacs/mdtypes/inputrec.h"
43 #include "gromacs/utility/fatalerror.h"
44 #include "gromacs/utility/smalloc.h"
45
46 gmx_enerdata_t::gmx_enerdata_t(int numEnergyGroups, int numFepLambdas) :
47     grpp(numEnergyGroups),
48     enerpart_lambda(numFepLambdas == 0 ? 0 : numFepLambdas + 1),
49     foreign_grpp(numEnergyGroups)
50 {
51 }
52
53 static real sum_v(int n, gmx::ArrayRef<const real> v)
54 {
55     real t;
56     int  i;
57
58     t = 0.0;
59     for (i = 0; (i < n); i++)
60     {
61         t = t + v[i];
62     }
63
64     return t;
65 }
66
67 void sum_epot(gmx_grppairener_t* grpp, real* epot)
68 {
69     int i;
70
71     /* Accumulate energies */
72     epot[F_COUL_SR] = sum_v(grpp->nener, grpp->ener[egCOULSR]);
73     epot[F_LJ]      = sum_v(grpp->nener, grpp->ener[egLJSR]);
74     epot[F_LJ14]    = sum_v(grpp->nener, grpp->ener[egLJ14]);
75     epot[F_COUL14]  = sum_v(grpp->nener, grpp->ener[egCOUL14]);
76
77     /* lattice part of LR doesnt belong to any group
78      * and has been added earlier
79      */
80     epot[F_BHAM] = sum_v(grpp->nener, grpp->ener[egBHAMSR]);
81
82     epot[F_EPOT] = 0;
83     for (i = 0; (i < F_EPOT); i++)
84     {
85         if (i != F_DISRESVIOL && i != F_ORIRESDEV)
86         {
87             epot[F_EPOT] += epot[i];
88         }
89     }
90 }
91
92 void sum_dhdl(gmx_enerdata_t* enerd, gmx::ArrayRef<const real> lambda, const t_lambda& fepvals)
93 {
94     int index;
95
96     enerd->dvdl_lin[efptVDW] += enerd->term[F_DVDL_VDW]; /* include dispersion correction */
97     enerd->term[F_DVDL] = 0.0;
98     for (int i = 0; i < efptNR; i++)
99     {
100         if (fepvals.separate_dvdl[i])
101         {
102             /* could this be done more readably/compactly? */
103             switch (i)
104             {
105                 case (efptMASS): index = F_DKDL; break;
106                 case (efptCOUL): index = F_DVDL_COUL; break;
107                 case (efptVDW): index = F_DVDL_VDW; break;
108                 case (efptBONDED): index = F_DVDL_BONDED; break;
109                 case (efptRESTRAINT): index = F_DVDL_RESTRAINT; break;
110                 default: index = F_DVDL; break;
111             }
112             enerd->term[index] = enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
113             if (debug)
114             {
115                 fprintf(debug, "dvdl-%s[%2d]: %f: non-linear %f + linear %f\n", efpt_names[i], i,
116                         enerd->term[index], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
117             }
118         }
119         else
120         {
121             enerd->term[F_DVDL] += enerd->dvdl_lin[i] + enerd->dvdl_nonlin[i];
122             if (debug)
123             {
124                 fprintf(debug, "dvd-%sl[%2d]: %f: non-linear %f + linear %f\n", efpt_names[0], i,
125                         enerd->term[F_DVDL], enerd->dvdl_nonlin[i], enerd->dvdl_lin[i]);
126             }
127         }
128     }
129
130     if (fepvals.separate_dvdl[efptBONDED])
131     {
132         enerd->term[F_DVDL_BONDED] += enerd->term[F_DVDL_CONSTR];
133     }
134     else
135     {
136         enerd->term[F_DVDL] += enerd->term[F_DVDL_CONSTR];
137     }
138
139     for (int i = 0; i < fepvals.n_lambda; i++)
140     {
141         /* note we are iterating over fepvals here!
142            For the current lam, dlam = 0 automatically,
143            so we don't need to add anything to the
144            enerd->enerpart_lambda[0] */
145
146         /* we don't need to worry about dvdl_lin contributions to dE at
147            current lambda, because the contributions to the current
148            lambda are automatically zeroed */
149
150         double& enerpart_lambda = enerd->enerpart_lambda[i + 1];
151
152         for (gmx::index j = 0; j < lambda.ssize(); j++)
153         {
154             /* Note that this loop is over all dhdl components, not just the separated ones */
155             const double dlam = fepvals.all_lambda[j][i] - lambda[j];
156
157             enerpart_lambda += dlam * enerd->dvdl_lin[j];
158
159             /* Constraints can not be evaluated at foreign lambdas, so we add
160              * a linear extrapolation. This is an approximation, but usually
161              * quite accurate since constraints change little between lambdas.
162              */
163             if ((j == efptBONDED && fepvals.separate_dvdl[efptBONDED])
164                 || (j == efptFEP && !fepvals.separate_dvdl[efptBONDED]))
165             {
166                 enerpart_lambda += dlam * enerd->term[F_DVDL_CONSTR];
167             }
168
169             if (j == efptMASS && !fepvals.separate_dvdl[j])
170             {
171                 enerpart_lambda += dlam * enerd->term[F_DKDL];
172             }
173
174             if (debug)
175             {
176                 fprintf(debug, "enerdiff lam %g: (%15s), non-linear %f linear %f*%f\n",
177                         fepvals.all_lambda[j][i], efpt_names[j],
178                         enerpart_lambda - enerd->enerpart_lambda[0], dlam, enerd->dvdl_lin[j]);
179             }
180         }
181     }
182
183     /* The constrain contribution is now included in other terms, so clear it */
184     enerd->term[F_DVDL_CONSTR] = 0;
185 }
186
187
188 void reset_foreign_enerdata(gmx_enerdata_t* enerd)
189 {
190     int i, j;
191
192     /* First reset all foreign energy components.  Foreign energies always called on
193        neighbor search steps */
194     for (i = 0; (i < egNR); i++)
195     {
196         for (j = 0; (j < enerd->grpp.nener); j++)
197         {
198             enerd->foreign_grpp.ener[i][j] = 0.0;
199         }
200     }
201
202     /* potential energy components */
203     for (i = 0; (i <= F_EPOT); i++)
204     {
205         enerd->foreign_term[i] = 0.0;
206     }
207 }
208
209 void reset_enerdata(gmx_enerdata_t* enerd)
210 {
211     int i, j;
212
213     /* First reset all energy components. */
214     for (i = 0; (i < egNR); i++)
215     {
216         for (j = 0; (j < enerd->grpp.nener); j++)
217         {
218             enerd->grpp.ener[i][j] = 0.0;
219         }
220     }
221     for (i = 0; i < efptNR; i++)
222     {
223         enerd->dvdl_lin[i]    = 0.0;
224         enerd->dvdl_nonlin[i] = 0.0;
225     }
226
227     /* Normal potential energy components */
228     for (i = 0; (i <= F_EPOT); i++)
229     {
230         enerd->term[i] = 0.0;
231     }
232     enerd->term[F_DVDL]           = 0.0;
233     enerd->term[F_DVDL_COUL]      = 0.0;
234     enerd->term[F_DVDL_VDW]       = 0.0;
235     enerd->term[F_DVDL_BONDED]    = 0.0;
236     enerd->term[F_DVDL_RESTRAINT] = 0.0;
237     enerd->term[F_DKDL]           = 0.0;
238     std::fill(enerd->enerpart_lambda.begin(), enerd->enerpart_lambda.end(), 0);
239     /* reset foreign energy data - separate function since we also call it elsewhere */
240     reset_foreign_enerdata(enerd);
241 }