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