Split lines with many copyright years
[alexxy/gromacs.git] / src / gromacs / mdlib / ebin.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) 2012,2014,2015,2017,2018 by the GROMACS development team.
7  * Copyright (c) 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 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
40
41 #include "ebin.h"
42
43 #include <cmath>
44 #include <cstring>
45
46 #include <algorithm>
47
48 #include "gromacs/math/units.h"
49 #include "gromacs/math/utilities.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/topology/ifunc.h"
52 #include "gromacs/trajectory/energyframe.h"
53 #include "gromacs/utility/basedefinitions.h"
54 #include "gromacs/utility/cstringutil.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/smalloc.h"
57 #include "gromacs/utility/stringutil.h"
58
59 t_ebin* mk_ebin()
60 {
61     t_ebin* eb;
62
63     snew(eb, 1);
64
65     return eb;
66 }
67
68 void done_ebin(t_ebin* eb)
69 {
70     for (int i = 0; i < eb->nener; i++)
71     {
72         sfree(eb->enm[i].name);
73         sfree(eb->enm[i].unit);
74     }
75     sfree(eb->e);
76     sfree(eb->e_sim);
77     sfree(eb->enm);
78     sfree(eb);
79 }
80
81 int get_ebin_space(t_ebin* eb, int nener, const char* const enm[], const char* unit)
82 {
83     int         index;
84     int         i, f;
85     const char* u;
86
87     index = eb->nener;
88     eb->nener += nener;
89     srenew(eb->e, eb->nener);
90     srenew(eb->e_sim, eb->nener);
91     srenew(eb->enm, eb->nener);
92     for (i = index; (i < eb->nener); i++)
93     {
94         eb->e[i].e        = 0;
95         eb->e[i].eav      = 0;
96         eb->e[i].esum     = 0;
97         eb->e_sim[i].e    = 0;
98         eb->e_sim[i].eav  = 0;
99         eb->e_sim[i].esum = 0;
100         eb->enm[i].name   = gmx_strdup(enm[i - index]);
101         if (unit != nullptr)
102         {
103             eb->enm[i].unit = gmx_strdup(unit);
104         }
105         else
106         {
107             /* Determine the unit from the longname.
108              * These units should have been defined in ifunc.c
109              * But even better would be if all interactions functions
110              * return energies and all non-interaction function
111              * entries would be removed from the ifunc array.
112              */
113             u = unit_energy;
114             for (f = 0; f < F_NRE; f++)
115             {
116                 if (strcmp(eb->enm[i].name, interaction_function[f].longname) == 0)
117                 {
118                     /* Only the terms in this list are not energies */
119                     switch (f)
120                     {
121                         case F_DISRESVIOL: u = unit_length; break;
122                         case F_ORIRESDEV: u = "obs"; break;
123                         case F_TEMP: u = unit_temp_K; break;
124                         case F_PDISPCORR:
125                         case F_PRES: u = unit_pres_bar; break;
126                     }
127                 }
128             }
129             eb->enm[i].unit = gmx_strdup(u);
130         }
131     }
132
133     return index;
134 }
135
136 // ICC 19 -O3 -msse2 generates wrong code. Lower optimization levels
137 // and other SIMD levels seem fine, however.
138 #if defined __ICC
139 #    pragma intel optimization_level 2
140 #endif
141 void add_ebin(t_ebin* eb, int entryIndex, int nener, const real ener[], gmx_bool bSum)
142 {
143     int       i, m;
144     double    e, invmm, diff;
145     t_energy *eg, *egs;
146
147     if ((entryIndex + nener > eb->nener) || (entryIndex < 0))
148     {
149         gmx_fatal(FARGS, "%s-%d: Energies out of range: entryIndex=%d nener=%d maxener=%d",
150                   __FILE__, __LINE__, entryIndex, nener, eb->nener);
151     }
152
153     eg = &(eb->e[entryIndex]);
154
155     for (i = 0; (i < nener); i++)
156     {
157         eg[i].e = ener[i];
158     }
159
160     if (bSum)
161     {
162         egs = &(eb->e_sim[entryIndex]);
163
164         m = eb->nsum;
165
166         if (m == 0)
167         {
168             for (i = 0; (i < nener); i++)
169             {
170                 eg[i].eav  = 0;
171                 eg[i].esum = ener[i];
172                 egs[i].esum += ener[i];
173             }
174         }
175         else
176         {
177             invmm = (1.0 / m) / (m + 1.0);
178
179             for (i = 0; (i < nener); i++)
180             {
181                 /* Value for this component */
182                 e = ener[i];
183
184                 /* first update sigma, then sum */
185                 diff = eg[i].esum - m * e;
186                 eg[i].eav += diff * diff * invmm;
187                 eg[i].esum += e;
188                 egs[i].esum += e;
189             }
190         }
191     }
192 }
193
194 // TODO It would be faster if this function was templated on both bSum
195 // and whether eb->nsum was zero, to lift the branches out of the loop
196 // over all possible energy terms, but that is true for a lot of the
197 // ebin and mdebin functionality, so we should do it all or nothing.
198 void add_ebin_indexed(t_ebin*                   eb,
199                       int                       entryIndex,
200                       gmx::ArrayRef<bool>       shouldUse,
201                       gmx::ArrayRef<const real> ener,
202                       gmx_bool                  bSum)
203 {
204
205     GMX_ASSERT(shouldUse.size() == ener.size(), "View sizes must match");
206     GMX_ASSERT(entryIndex + std::count(shouldUse.begin(), shouldUse.end(), true) <= eb->nener,
207                gmx::formatString("Energies out of range: entryIndex=%d nener=%td maxener=%d", entryIndex,
208                                  std::count(shouldUse.begin(), shouldUse.end(), true), eb->nener)
209                        .c_str());
210     GMX_ASSERT(entryIndex >= 0, "Must have non-negative entry");
211
212     const int    m              = eb->nsum;
213     const double invmm          = (m == 0) ? 0 : (1.0 / m) / (m + 1.0);
214     t_energy*    energyEntry    = &(eb->e[entryIndex]);
215     t_energy*    simEnergyEntry = &(eb->e_sim[entryIndex]);
216     auto         shouldUseIter  = shouldUse.begin();
217     for (const auto& theEnergy : ener)
218     {
219         if (*shouldUseIter)
220         {
221             energyEntry->e = theEnergy;
222             if (bSum)
223             {
224                 if (m == 0)
225                 {
226                     energyEntry->eav  = 0;
227                     energyEntry->esum = theEnergy;
228                     simEnergyEntry->esum += theEnergy;
229                 }
230                 else
231                 {
232                     /* first update sigma, then sum */
233                     double diff = energyEntry->esum - m * theEnergy;
234                     energyEntry->eav += diff * diff * invmm;
235                     energyEntry->esum += theEnergy;
236                     simEnergyEntry->esum += theEnergy;
237                 }
238                 ++simEnergyEntry;
239             }
240             ++energyEntry;
241         }
242         ++shouldUseIter;
243     }
244 }
245
246 void ebin_increase_count(int increment, t_ebin* eb, gmx_bool bSum)
247 {
248     eb->nsteps += increment;
249     eb->nsteps_sim += increment;
250
251     if (bSum)
252     {
253         eb->nsum += increment;
254         eb->nsum_sim += increment;
255     }
256 }
257
258 void reset_ebin_sums(t_ebin* eb)
259 {
260     eb->nsteps = 0;
261     eb->nsum   = 0;
262     /* The actual sums are cleared when the next frame is stored */
263 }
264
265 void pr_ebin(FILE* fp, t_ebin* eb, int entryIndex, int nener, int nperline, int prmode, gmx_bool bPrHead)
266 {
267     int  i, j, i0;
268     int  rc;
269     char buf[30];
270
271     rc = 0;
272
273     if (entryIndex < 0 || entryIndex > eb->nener)
274     {
275         gmx_fatal(FARGS, "Invalid entryIndex in pr_ebin: %d", entryIndex);
276     }
277     int start = entryIndex;
278     if (nener > eb->nener)
279     {
280         gmx_fatal(FARGS, "Invalid nener in pr_ebin: %d", nener);
281     }
282     int end = eb->nener;
283     if (nener != -1)
284     {
285         end = entryIndex + nener;
286     }
287     for (i = start; (i < end) && rc >= 0;)
288     {
289         if (bPrHead)
290         {
291             i0 = i;
292             for (j = 0; (j < nperline) && (i < end) && rc >= 0; j++, i++)
293             {
294                 if (strncmp(eb->enm[i].name, "Pres", 4) == 0)
295                 {
296                     /* Print the pressure unit to avoid confusion */
297                     sprintf(buf, "%s (%s)", eb->enm[i].name, unit_pres_bar);
298                     rc = fprintf(fp, "%15s", buf);
299                 }
300                 else
301                 {
302                     rc = fprintf(fp, "%15s", eb->enm[i].name);
303                 }
304             }
305
306             if (rc >= 0)
307             {
308                 rc = fprintf(fp, "\n");
309             }
310
311             i = i0;
312         }
313         for (j = 0; (j < nperline) && (i < end) && rc >= 0; j++, i++)
314         {
315             switch (prmode)
316             {
317                 case eprNORMAL: rc = fprintf(fp, "   %12.5e", eb->e[i].e); break;
318                 case eprAVER:
319                     if (eb->nsum_sim > 0)
320                     {
321                         rc = fprintf(fp, "   %12.5e", eb->e_sim[i].esum / eb->nsum_sim);
322                     }
323                     else
324                     {
325                         rc = fprintf(fp, "    %-12s", "N/A");
326                     }
327                     break;
328                 default: gmx_fatal(FARGS, "Invalid print mode %d in pr_ebin", prmode);
329             }
330         }
331         if (rc >= 0)
332         {
333             rc = fprintf(fp, "\n");
334         }
335     }
336     if (rc < 0)
337     {
338         gmx_fatal(FARGS, "Cannot write to logfile; maybe you are out of disk space?");
339     }
340 }