Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / mdlib / stat.c
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, 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 <stdio.h>
40 #include <string.h>
41
42 #include "gromacs/fileio/gmxfio.h"
43 #include "gromacs/fileio/trnio.h"
44 #include "gromacs/fileio/xtcio.h"
45 #include "gromacs/legacyheaders/checkpoint.h"
46 #include "gromacs/legacyheaders/constr.h"
47 #include "gromacs/legacyheaders/domdec.h"
48 #include "gromacs/legacyheaders/force.h"
49 #include "gromacs/legacyheaders/md_support.h"
50 #include "gromacs/legacyheaders/mdrun.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gromacs/legacyheaders/network.h"
53 #include "gromacs/legacyheaders/rbin.h"
54 #include "gromacs/legacyheaders/sim_util.h"
55 #include "gromacs/legacyheaders/tgroup.h"
56 #include "gromacs/legacyheaders/txtdump.h"
57 #include "gromacs/legacyheaders/typedefs.h"
58 #include "gromacs/legacyheaders/vcm.h"
59 #include "gromacs/legacyheaders/types/commrec.h"
60 #include "gromacs/math/utilities.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/futil.h"
64 #include "gromacs/utility/smalloc.h"
65
66 typedef struct gmx_global_stat
67 {
68     t_bin *rb;
69     int   *itc0;
70     int   *itc1;
71 } t_gmx_global_stat;
72
73 gmx_global_stat_t global_stat_init(t_inputrec *ir)
74 {
75     gmx_global_stat_t gs;
76
77     snew(gs, 1);
78
79     gs->rb = mk_bin();
80     snew(gs->itc0, ir->opts.ngtc);
81     snew(gs->itc1, ir->opts.ngtc);
82
83     return gs;
84 }
85
86 void global_stat_destroy(gmx_global_stat_t gs)
87 {
88     destroy_bin(gs->rb);
89     sfree(gs->itc0);
90     sfree(gs->itc1);
91     sfree(gs);
92 }
93
94 static int filter_enerdterm(real *afrom, gmx_bool bToBuffer, real *ato,
95                             gmx_bool bTemp, gmx_bool bPres, gmx_bool bEner)
96 {
97     int i, to, from;
98
99     from = 0;
100     to   = 0;
101     for (i = 0; i < F_NRE; i++)
102     {
103         if (bToBuffer)
104         {
105             from = i;
106         }
107         else
108         {
109             to = i;
110         }
111         switch (i)
112         {
113             case F_EKIN:
114             case F_TEMP:
115             case F_DKDL:
116                 if (bTemp)
117                 {
118                     ato[to++] = afrom[from++];
119                 }
120                 break;
121             case F_PRES:
122             case F_PDISPCORR:
123                 if (bPres)
124                 {
125                     ato[to++] = afrom[from++];
126                 }
127                 break;
128             default:
129                 if (bEner)
130                 {
131                     ato[to++] = afrom[from++];
132                 }
133                 break;
134         }
135     }
136
137     return to;
138 }
139
140 void global_stat(FILE *fplog, gmx_global_stat_t gs,
141                  t_commrec *cr, gmx_enerdata_t *enerd,
142                  tensor fvir, tensor svir, rvec mu_tot,
143                  t_inputrec *inputrec,
144                  gmx_ekindata_t *ekind, gmx_constr_t constr,
145                  t_vcm *vcm,
146                  int nsig, real *sig,
147                  gmx_mtop_t *top_global, t_state *state_local,
148                  gmx_bool bSumEkinhOld, int flags)
149 /* instead of current system, gmx_booleans for summing virial, kinetic energy, and other terms */
150 {
151     t_bin     *rb;
152     int       *itc0, *itc1;
153     int        ie    = 0, ifv = 0, isv = 0, irmsd = 0, imu = 0;
154     int        idedl = 0, idvdll = 0, idvdlnl = 0, iepl = 0, icm = 0, imass = 0, ica = 0, inb = 0;
155     int        isig  = -1;
156     int        icj   = -1, ici = -1, icx = -1;
157     int        inn[egNR];
158     real       copyenerd[F_NRE];
159     int        nener, j;
160     real      *rmsd_data = NULL;
161     double     nb;
162     gmx_bool   bVV, bTemp, bEner, bPres, bConstrVir, bEkinAveVel, bFirstIterate, bReadEkin;
163
164     bVV           = EI_VV(inputrec->eI);
165     bTemp         = flags & CGLO_TEMPERATURE;
166     bEner         = flags & CGLO_ENERGY;
167     bPres         = (flags & CGLO_PRESSURE);
168     bConstrVir    = (flags & CGLO_CONSTRAINT);
169     bFirstIterate = (flags & CGLO_FIRSTITERATE);
170     bEkinAveVel   = (inputrec->eI == eiVV || (inputrec->eI == eiVVAK && bPres));
171     bReadEkin     = (flags & CGLO_READEKIN);
172
173     rb   = gs->rb;
174     itc0 = gs->itc0;
175     itc1 = gs->itc1;
176
177
178     reset_bin(rb);
179     /* This routine copies all the data to be summed to one big buffer
180      * using the t_bin struct.
181      */
182
183     /* First, we neeed to identify which enerd->term should be
184        communicated.  Temperature and pressure terms should only be
185        communicated and summed when they need to be, to avoid repeating
186        the sums and overcounting. */
187
188     nener = filter_enerdterm(enerd->term, TRUE, copyenerd, bTemp, bPres, bEner);
189
190     /* First, the data that needs to be communicated with velocity verlet every time
191        This is just the constraint virial.*/
192     if (bConstrVir)
193     {
194         isv = add_binr(rb, DIM*DIM, svir[0]);
195         where();
196     }
197
198 /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
199     if (bTemp || !bVV)
200     {
201         if (ekind)
202         {
203             for (j = 0; (j < inputrec->opts.ngtc); j++)
204             {
205                 if (bSumEkinhOld)
206                 {
207                     itc0[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
208                 }
209                 if (bEkinAveVel && !bReadEkin)
210                 {
211                     itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinf[0]);
212                 }
213                 else if (!bReadEkin)
214                 {
215                     itc1[j] = add_binr(rb, DIM*DIM, ekind->tcstat[j].ekinh[0]);
216                 }
217             }
218             /* these probably need to be put into one of these categories */
219             where();
220             idedl = add_binr(rb, 1, &(ekind->dekindl));
221             where();
222             ica   = add_binr(rb, 1, &(ekind->cosacc.mvcos));
223             where();
224         }
225     }
226     where();
227
228     if ((bPres || !bVV) && bFirstIterate)
229     {
230         ifv = add_binr(rb, DIM*DIM, fvir[0]);
231     }
232
233
234     if (bEner)
235     {
236         where();
237         if (bFirstIterate)
238         {
239             ie  = add_binr(rb, nener, copyenerd);
240         }
241         where();
242         if (constr)
243         {
244             rmsd_data = constr_rmsd_data(constr);
245             if (rmsd_data)
246             {
247                 irmsd = add_binr(rb, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
248             }
249         }
250         if (!NEED_MUTOT(*inputrec))
251         {
252             imu = add_binr(rb, DIM, mu_tot);
253             where();
254         }
255
256         if (bFirstIterate)
257         {
258             for (j = 0; (j < egNR); j++)
259             {
260                 inn[j] = add_binr(rb, enerd->grpp.nener, enerd->grpp.ener[j]);
261             }
262             where();
263             if (inputrec->efep != efepNO)
264             {
265                 idvdll  = add_bind(rb, efptNR, enerd->dvdl_lin);
266                 idvdlnl = add_bind(rb, efptNR, enerd->dvdl_nonlin);
267                 if (enerd->n_lambda > 0)
268                 {
269                     iepl = add_bind(rb, enerd->n_lambda, enerd->enerpart_lambda);
270                 }
271             }
272         }
273     }
274
275     if (vcm)
276     {
277         icm   = add_binr(rb, DIM*vcm->nr, vcm->group_p[0]);
278         where();
279         imass = add_binr(rb, vcm->nr, vcm->group_mass);
280         where();
281         if (vcm->mode == ecmANGULAR)
282         {
283             icj   = add_binr(rb, DIM*vcm->nr, vcm->group_j[0]);
284             where();
285             icx   = add_binr(rb, DIM*vcm->nr, vcm->group_x[0]);
286             where();
287             ici   = add_binr(rb, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
288             where();
289         }
290     }
291
292     if (DOMAINDECOMP(cr))
293     {
294         nb  = cr->dd->nbonded_local;
295         inb = add_bind(rb, 1, &nb);
296     }
297     where();
298     if (nsig > 0)
299     {
300         isig = add_binr(rb, nsig, sig);
301     }
302
303     /* Global sum it all */
304     if (debug)
305     {
306         fprintf(debug, "Summing %d energies\n", rb->maxreal);
307     }
308     sum_bin(rb, cr);
309     where();
310
311     /* Extract all the data locally */
312
313     if (bConstrVir)
314     {
315         extract_binr(rb, isv, DIM*DIM, svir[0]);
316     }
317
318     /* We need the force virial and the kinetic energy for the first time through with velocity verlet */
319     if (bTemp || !bVV)
320     {
321         if (ekind)
322         {
323             for (j = 0; (j < inputrec->opts.ngtc); j++)
324             {
325                 if (bSumEkinhOld)
326                 {
327                     extract_binr(rb, itc0[j], DIM*DIM, ekind->tcstat[j].ekinh_old[0]);
328                 }
329                 if (bEkinAveVel && !bReadEkin)
330                 {
331                     extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinf[0]);
332                 }
333                 else if (!bReadEkin)
334                 {
335                     extract_binr(rb, itc1[j], DIM*DIM, ekind->tcstat[j].ekinh[0]);
336                 }
337             }
338             extract_binr(rb, idedl, 1, &(ekind->dekindl));
339             extract_binr(rb, ica, 1, &(ekind->cosacc.mvcos));
340             where();
341         }
342     }
343     if ((bPres || !bVV) && bFirstIterate)
344     {
345         extract_binr(rb, ifv, DIM*DIM, fvir[0]);
346     }
347
348     if (bEner)
349     {
350         if (bFirstIterate)
351         {
352             extract_binr(rb, ie, nener, copyenerd);
353             if (rmsd_data)
354             {
355                 extract_binr(rb, irmsd, inputrec->eI == eiSD2 ? 3 : 2, rmsd_data);
356             }
357             if (!NEED_MUTOT(*inputrec))
358             {
359                 extract_binr(rb, imu, DIM, mu_tot);
360             }
361
362             for (j = 0; (j < egNR); j++)
363             {
364                 extract_binr(rb, inn[j], enerd->grpp.nener, enerd->grpp.ener[j]);
365             }
366             if (inputrec->efep != efepNO)
367             {
368                 extract_bind(rb, idvdll, efptNR, enerd->dvdl_lin);
369                 extract_bind(rb, idvdlnl, efptNR, enerd->dvdl_nonlin);
370                 if (enerd->n_lambda > 0)
371                 {
372                     extract_bind(rb, iepl, enerd->n_lambda, enerd->enerpart_lambda);
373                 }
374             }
375             if (DOMAINDECOMP(cr))
376             {
377                 extract_bind(rb, inb, 1, &nb);
378                 if ((int)(nb + 0.5) != cr->dd->nbonded_global)
379                 {
380                     dd_print_missing_interactions(fplog, cr, (int)(nb + 0.5), top_global, state_local);
381                 }
382             }
383             where();
384
385             filter_enerdterm(copyenerd, FALSE, enerd->term, bTemp, bPres, bEner);
386         }
387     }
388
389     if (vcm)
390     {
391         extract_binr(rb, icm, DIM*vcm->nr, vcm->group_p[0]);
392         where();
393         extract_binr(rb, imass, vcm->nr, vcm->group_mass);
394         where();
395         if (vcm->mode == ecmANGULAR)
396         {
397             extract_binr(rb, icj, DIM*vcm->nr, vcm->group_j[0]);
398             where();
399             extract_binr(rb, icx, DIM*vcm->nr, vcm->group_x[0]);
400             where();
401             extract_binr(rb, ici, DIM*DIM*vcm->nr, vcm->group_i[0][0]);
402             where();
403         }
404     }
405
406     if (nsig > 0)
407     {
408         extract_binr(rb, isig, nsig, sig);
409     }
410     where();
411 }
412
413 int do_per_step(gmx_int64_t step, gmx_int64_t nstep)
414 {
415     if (nstep != 0)
416     {
417         return ((step % nstep) == 0);
418     }
419     else
420     {
421         return 0;
422     }
423 }