Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / mdebin.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 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <float.h>
42 #include <stdlib.h>
43 #include <string.h>
44
45 #include "typedefs.h"
46 #include "mdebin.h"
47 #include "gromacs/math/units.h"
48 #include "gromacs/fileio/enxio.h"
49 #include "gromacs/math/vec.h"
50 #include "disre.h"
51 #include "network.h"
52 #include "names.h"
53 #include "orires.h"
54 #include "constr.h"
55 #include "gromacs/topology/mtop_util.h"
56 #include "gromacs/fileio/xvgr.h"
57 #include "gromacs/fileio/gmxfio.h"
58 #include "macros.h"
59 #include "mdrun.h"
60 #include "mdebin_bar.h"
61
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/utility/smalloc.h"
64
65 static const char *conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
66
67 static const char *boxs_nm[] = { "Box-X", "Box-Y", "Box-Z" };
68
69 static const char *tricl_boxs_nm[] = {
70     "Box-XX", "Box-YY", "Box-ZZ",
71     "Box-YX", "Box-ZX", "Box-ZY"
72 };
73
74 static const char *vol_nm[] = { "Volume" };
75
76 static const char *dens_nm[] = {"Density" };
77
78 static const char *pv_nm[] = {"pV" };
79
80 static const char *enthalpy_nm[] = {"Enthalpy" };
81
82 static const char *boxvel_nm[] = {
83     "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
84     "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
85 };
86
87 #define NBOXS asize(boxs_nm)
88 #define NTRICLBOXS asize(tricl_boxs_nm)
89
90 t_mdebin *init_mdebin(ener_file_t       fp_ene,
91                       const gmx_mtop_t *mtop,
92                       const t_inputrec *ir,
93                       FILE             *fp_dhdl)
94 {
95     const char         *ener_nm[F_NRE];
96     static const char  *vir_nm[] = {
97         "Vir-XX", "Vir-XY", "Vir-XZ",
98         "Vir-YX", "Vir-YY", "Vir-YZ",
99         "Vir-ZX", "Vir-ZY", "Vir-ZZ"
100     };
101     static const char  *sv_nm[] = {
102         "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
103         "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
104         "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ"
105     };
106     static const char  *fv_nm[] = {
107         "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
108         "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
109         "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ"
110     };
111     static const char  *pres_nm[] = {
112         "Pres-XX", "Pres-XY", "Pres-XZ",
113         "Pres-YX", "Pres-YY", "Pres-YZ",
114         "Pres-ZX", "Pres-ZY", "Pres-ZZ"
115     };
116     static const char  *surft_nm[] = {
117         "#Surf*SurfTen"
118     };
119     static const char  *mu_nm[] = {
120         "Mu-X", "Mu-Y", "Mu-Z"
121     };
122     static const char  *vcos_nm[] = {
123         "2CosZ*Vel-X"
124     };
125     static const char  *visc_nm[] = {
126         "1/Viscosity"
127     };
128     static const char  *baro_nm[] = {
129         "Barostat"
130     };
131
132     char              **grpnms;
133     const gmx_groups_t *groups;
134     char              **gnm;
135     char                buf[256];
136     const char         *bufi;
137     t_mdebin           *md;
138     int                 i, j, ni, nj, n, nh, k, kk, ncon, nset;
139     gmx_bool            bBHAM, bNoseHoover, b14;
140
141     snew(md, 1);
142
143     md->bVir   = TRUE;
144     md->bPress = TRUE;
145     md->bSurft = TRUE;
146     md->bMu    = TRUE;
147
148     if (EI_DYNAMICS(ir->eI))
149     {
150         md->delta_t = ir->delta_t;
151     }
152     else
153     {
154         md->delta_t = 0;
155     }
156
157     groups = &mtop->groups;
158
159     bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
160     b14   = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
161              gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
162
163     ncon           = gmx_mtop_ftype_count(mtop, F_CONSTR);
164     nset           = gmx_mtop_ftype_count(mtop, F_SETTLE);
165     md->bConstr    = (ncon > 0 || nset > 0);
166     md->bConstrVir = FALSE;
167     if (md->bConstr)
168     {
169         if (ncon > 0 && ir->eConstrAlg == econtLINCS)
170         {
171             if (ir->eI == eiSD2)
172             {
173                 md->nCrmsd = 2;
174             }
175             else
176             {
177                 md->nCrmsd = 1;
178             }
179         }
180         md->bConstrVir = (getenv("GMX_CONSTRAINTVIR") != NULL);
181     }
182     else
183     {
184         md->nCrmsd = 0;
185     }
186
187     /* Energy monitoring */
188     for (i = 0; i < egNR; i++)
189     {
190         md->bEInd[i] = FALSE;
191     }
192
193     for (i = 0; i < F_NRE; i++)
194     {
195         md->bEner[i] = FALSE;
196         if (i == F_LJ)
197         {
198             md->bEner[i] = !bBHAM;
199         }
200         else if (i == F_BHAM)
201         {
202             md->bEner[i] = bBHAM;
203         }
204         else if (i == F_EQM)
205         {
206             md->bEner[i] = ir->bQMMM;
207         }
208         else if (i == F_COUL_LR)
209         {
210             md->bEner[i] = (ir->rcoulomb > ir->rlist);
211         }
212         else if (i == F_LJ_LR)
213         {
214             md->bEner[i] = (!bBHAM && ir->rvdw > ir->rlist);
215         }
216         else if (i == F_BHAM_LR)
217         {
218             md->bEner[i] = (bBHAM && ir->rvdw > ir->rlist);
219         }
220         else if (i == F_RF_EXCL)
221         {
222             md->bEner[i] = (EEL_RF(ir->coulombtype) && ir->coulombtype != eelRF_NEC && ir->cutoff_scheme == ecutsGROUP);
223         }
224         else if (i == F_COUL_RECIP)
225         {
226             md->bEner[i] = EEL_FULL(ir->coulombtype);
227         }
228         else if (i == F_LJ_RECIP)
229         {
230             md->bEner[i] = EVDW_PME(ir->vdwtype);
231         }
232         else if (i == F_LJ14)
233         {
234             md->bEner[i] = b14;
235         }
236         else if (i == F_COUL14)
237         {
238             md->bEner[i] = b14;
239         }
240         else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
241         {
242             md->bEner[i] = FALSE;
243         }
244         else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
245                  (i == F_DVDL_VDW  && ir->fepvals->separate_dvdl[efptVDW]) ||
246                  (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
247                  (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
248                  (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
249                  (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
250         {
251             md->bEner[i] = (ir->efep != efepNO);
252         }
253         else if ((interaction_function[i].flags & IF_VSITE) ||
254                  (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
255         {
256             md->bEner[i] = FALSE;
257         }
258         else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES)  || (i == F_EQM))
259         {
260             md->bEner[i] = TRUE;
261         }
262         else if ((i == F_GBPOL) && ir->implicit_solvent == eisGBSA)
263         {
264             md->bEner[i] = TRUE;
265         }
266         else if ((i == F_NPSOLVATION) && ir->implicit_solvent == eisGBSA && (ir->sa_algorithm != esaNO))
267         {
268             md->bEner[i] = TRUE;
269         }
270         else if ((i == F_GB12) || (i == F_GB13) || (i == F_GB14))
271         {
272             md->bEner[i] = FALSE;
273         }
274         else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
275         {
276             md->bEner[i] = EI_DYNAMICS(ir->eI);
277         }
278         else if (i == F_DISPCORR || i == F_PDISPCORR)
279         {
280             md->bEner[i] = (ir->eDispCorr != edispcNO);
281         }
282         else if (i == F_DISRESVIOL)
283         {
284             md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
285         }
286         else if (i == F_ORIRESDEV)
287         {
288             md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
289         }
290         else if (i == F_CONNBONDS)
291         {
292             md->bEner[i] = FALSE;
293         }
294         else if (i == F_COM_PULL)
295         {
296             md->bEner[i] = (ir->ePull == epullUMBRELLA || ir->ePull == epullCONST_F || ir->bRot);
297         }
298         else if (i == F_ECONSERVED)
299         {
300             md->bEner[i] = ((ir->etc == etcNOSEHOOVER || ir->etc == etcVRESCALE) &&
301                             (ir->epc == epcNO || ir->epc == epcMTTK));
302         }
303         else
304         {
305             md->bEner[i] = (gmx_mtop_ftype_count(mtop, i) > 0);
306         }
307     }
308
309     /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
310     if (ir->bAdress && !debug)
311     {
312         for (i = 0; i < F_NRE; i++)
313         {
314             md->bEner[i] = FALSE;
315             if (i == F_EKIN)
316             {
317                 md->bEner[i] = TRUE;
318             }
319             if (i == F_TEMP)
320             {
321                 md->bEner[i] = TRUE;
322             }
323         }
324         md->bVir   = FALSE;
325         md->bPress = FALSE;
326         md->bSurft = FALSE;
327         md->bMu    = FALSE;
328     }
329
330     md->f_nre = 0;
331     for (i = 0; i < F_NRE; i++)
332     {
333         if (md->bEner[i])
334         {
335             ener_nm[md->f_nre] = interaction_function[i].longname;
336             md->f_nre++;
337         }
338     }
339
340     md->epc            = ir->epc;
341     md->bDiagPres      = !TRICLINIC(ir->ref_p);
342     md->ref_p          = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
343     md->bTricl         = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
344     md->bDynBox        = DYNAMIC_BOX(*ir);
345     md->etc            = ir->etc;
346     md->bNHC_trotter   = IR_NVT_TROTTER(ir);
347     md->bPrintNHChains = ir->bPrintNHChains;
348     md->bMTTK          = (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir));
349     md->bMu            = NEED_MUTOT(*ir);
350
351     md->ebin  = mk_ebin();
352     /* Pass NULL for unit to let get_ebin_space determine the units
353      * for interaction_function[i].longname
354      */
355     md->ie    = get_ebin_space(md->ebin, md->f_nre, ener_nm, NULL);
356     if (md->nCrmsd)
357     {
358         /* This should be called directly after the call for md->ie,
359          * such that md->iconrmsd follows directly in the list.
360          */
361         md->iconrmsd = get_ebin_space(md->ebin, md->nCrmsd, conrmsd_nm, "");
362     }
363     if (md->bDynBox)
364     {
365         md->ib    = get_ebin_space(md->ebin,
366                                    md->bTricl ? NTRICLBOXS : NBOXS,
367                                    md->bTricl ? tricl_boxs_nm : boxs_nm,
368                                    unit_length);
369         md->ivol  = get_ebin_space(md->ebin, 1, vol_nm,  unit_volume);
370         md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
371         if (md->bDiagPres)
372         {
373             md->ipv       = get_ebin_space(md->ebin, 1, pv_nm,   unit_energy);
374             md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm,   unit_energy);
375         }
376     }
377     if (md->bConstrVir)
378     {
379         md->isvir = get_ebin_space(md->ebin, asize(sv_nm), sv_nm, unit_energy);
380         md->ifvir = get_ebin_space(md->ebin, asize(fv_nm), fv_nm, unit_energy);
381     }
382     if (md->bVir)
383     {
384         md->ivir   = get_ebin_space(md->ebin, asize(vir_nm), vir_nm, unit_energy);
385     }
386     if (md->bPress)
387     {
388         md->ipres  = get_ebin_space(md->ebin, asize(pres_nm), pres_nm, unit_pres_bar);
389     }
390     if (md->bSurft)
391     {
392         md->isurft = get_ebin_space(md->ebin, asize(surft_nm), surft_nm,
393                                     unit_surft_bar);
394     }
395     if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
396     {
397         md->ipc = get_ebin_space(md->ebin, md->bTricl ? 6 : 3,
398                                  boxvel_nm, unit_vel);
399     }
400     if (md->bMu)
401     {
402         md->imu    = get_ebin_space(md->ebin, asize(mu_nm), mu_nm, unit_dipole_D);
403     }
404     if (ir->cos_accel != 0)
405     {
406         md->ivcos = get_ebin_space(md->ebin, asize(vcos_nm), vcos_nm, unit_vel);
407         md->ivisc = get_ebin_space(md->ebin, asize(visc_nm), visc_nm,
408                                    unit_invvisc_SI);
409     }
410
411     /* Energy monitoring */
412     for (i = 0; i < egNR; i++)
413     {
414         md->bEInd[i] = FALSE;
415     }
416     md->bEInd[egCOULSR] = TRUE;
417     md->bEInd[egLJSR  ] = TRUE;
418
419     if (ir->rcoulomb > ir->rlist)
420     {
421         md->bEInd[egCOULLR] = TRUE;
422     }
423     if (!bBHAM)
424     {
425         if (ir->rvdw > ir->rlist)
426         {
427             md->bEInd[egLJLR]   = TRUE;
428         }
429     }
430     else
431     {
432         md->bEInd[egLJSR]   = FALSE;
433         md->bEInd[egBHAMSR] = TRUE;
434         if (ir->rvdw > ir->rlist)
435         {
436             md->bEInd[egBHAMLR]   = TRUE;
437         }
438     }
439     if (b14)
440     {
441         md->bEInd[egLJ14]   = TRUE;
442         md->bEInd[egCOUL14] = TRUE;
443     }
444     md->nEc = 0;
445     for (i = 0; (i < egNR); i++)
446     {
447         if (md->bEInd[i])
448         {
449             md->nEc++;
450         }
451     }
452
453     n = groups->grps[egcENER].nr;
454     /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
455     if (!ir->bAdress)
456     {
457         /*standard simulation*/
458         md->nEg = n;
459         md->nE  = (n*(n+1))/2;
460     }
461     else if (!debug)
462     {
463         /*AdResS simulation*/
464         md->nU    = 0;
465         md->nEg   = 0;
466         md->nE    = 0;
467         md->nEc   = 0;
468         md->isvir = FALSE;
469     }
470     snew(md->igrp, md->nE);
471     if (md->nE > 1)
472     {
473         n = 0;
474         snew(gnm, md->nEc);
475         for (k = 0; (k < md->nEc); k++)
476         {
477             snew(gnm[k], STRLEN);
478         }
479         for (i = 0; (i < groups->grps[egcENER].nr); i++)
480         {
481             ni = groups->grps[egcENER].nm_ind[i];
482             for (j = i; (j < groups->grps[egcENER].nr); j++)
483             {
484                 nj = groups->grps[egcENER].nm_ind[j];
485                 for (k = kk = 0; (k < egNR); k++)
486                 {
487                     if (md->bEInd[k])
488                     {
489                         sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
490                                 *(groups->grpname[ni]), *(groups->grpname[nj]));
491                         kk++;
492                     }
493                 }
494                 md->igrp[n] = get_ebin_space(md->ebin, md->nEc,
495                                              (const char **)gnm, unit_energy);
496                 n++;
497             }
498         }
499         for (k = 0; (k < md->nEc); k++)
500         {
501             sfree(gnm[k]);
502         }
503         sfree(gnm);
504
505         if (n != md->nE)
506         {
507             gmx_incons("Number of energy terms wrong");
508         }
509     }
510
511     md->nTC  = groups->grps[egcTC].nr;
512     md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
513     if (md->bMTTK)
514     {
515         md->nTCP = 1;  /* assume only one possible coupling system for barostat
516                           for now */
517     }
518     else
519     {
520         md->nTCP = 0;
521     }
522     if (md->etc == etcNOSEHOOVER)
523     {
524         if (md->bNHC_trotter)
525         {
526             md->mde_n = 2*md->nNHC*md->nTC;
527         }
528         else
529         {
530             md->mde_n = 2*md->nTC;
531         }
532         if (md->epc == epcMTTK)
533         {
534             md->mdeb_n = 2*md->nNHC*md->nTCP;
535         }
536     }
537     else
538     {
539         md->mde_n  = md->nTC;
540         md->mdeb_n = 0;
541     }
542
543     snew(md->tmp_r, md->mde_n);
544     snew(md->tmp_v, md->mde_n);
545     snew(md->grpnms, md->mde_n);
546     grpnms = md->grpnms;
547
548     for (i = 0; (i < md->nTC); i++)
549     {
550         ni = groups->grps[egcTC].nm_ind[i];
551         sprintf(buf, "T-%s", *(groups->grpname[ni]));
552         grpnms[i] = strdup(buf);
553     }
554     md->itemp = get_ebin_space(md->ebin, md->nTC, (const char **)grpnms,
555                                unit_temp_K);
556
557     if (md->etc == etcNOSEHOOVER)
558     {
559         if (md->bPrintNHChains)
560         {
561             if (md->bNHC_trotter)
562             {
563                 for (i = 0; (i < md->nTC); i++)
564                 {
565                     ni   = groups->grps[egcTC].nm_ind[i];
566                     bufi = *(groups->grpname[ni]);
567                     for (j = 0; (j < md->nNHC); j++)
568                     {
569                         sprintf(buf, "Xi-%d-%s", j, bufi);
570                         grpnms[2*(i*md->nNHC+j)] = strdup(buf);
571                         sprintf(buf, "vXi-%d-%s", j, bufi);
572                         grpnms[2*(i*md->nNHC+j)+1] = strdup(buf);
573                     }
574                 }
575                 md->itc = get_ebin_space(md->ebin, md->mde_n,
576                                          (const char **)grpnms, unit_invtime);
577                 if (md->bMTTK)
578                 {
579                     for (i = 0; (i < md->nTCP); i++)
580                     {
581                         bufi = baro_nm[0];  /* All barostat DOF's together for now. */
582                         for (j = 0; (j < md->nNHC); j++)
583                         {
584                             sprintf(buf, "Xi-%d-%s", j, bufi);
585                             grpnms[2*(i*md->nNHC+j)] = strdup(buf);
586                             sprintf(buf, "vXi-%d-%s", j, bufi);
587                             grpnms[2*(i*md->nNHC+j)+1] = strdup(buf);
588                         }
589                     }
590                     md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
591                                               (const char **)grpnms, unit_invtime);
592                 }
593             }
594             else
595             {
596                 for (i = 0; (i < md->nTC); i++)
597                 {
598                     ni   = groups->grps[egcTC].nm_ind[i];
599                     bufi = *(groups->grpname[ni]);
600                     sprintf(buf, "Xi-%s", bufi);
601                     grpnms[2*i] = strdup(buf);
602                     sprintf(buf, "vXi-%s", bufi);
603                     grpnms[2*i+1] = strdup(buf);
604                 }
605                 md->itc = get_ebin_space(md->ebin, md->mde_n,
606                                          (const char **)grpnms, unit_invtime);
607             }
608         }
609     }
610     else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
611              md->etc == etcVRESCALE)
612     {
613         for (i = 0; (i < md->nTC); i++)
614         {
615             ni = groups->grps[egcTC].nm_ind[i];
616             sprintf(buf, "Lamb-%s", *(groups->grpname[ni]));
617             grpnms[i] = strdup(buf);
618         }
619         md->itc = get_ebin_space(md->ebin, md->mde_n, (const char **)grpnms, "");
620     }
621
622     sfree(grpnms);
623
624
625     md->nU = groups->grps[egcACC].nr;
626     if (md->nU > 1)
627     {
628         snew(grpnms, 3*md->nU);
629         for (i = 0; (i < md->nU); i++)
630         {
631             ni = groups->grps[egcACC].nm_ind[i];
632             sprintf(buf, "Ux-%s", *(groups->grpname[ni]));
633             grpnms[3*i+XX] = strdup(buf);
634             sprintf(buf, "Uy-%s", *(groups->grpname[ni]));
635             grpnms[3*i+YY] = strdup(buf);
636             sprintf(buf, "Uz-%s", *(groups->grpname[ni]));
637             grpnms[3*i+ZZ] = strdup(buf);
638         }
639         md->iu = get_ebin_space(md->ebin, 3*md->nU, (const char **)grpnms, unit_vel);
640         sfree(grpnms);
641     }
642
643     if (fp_ene)
644     {
645         do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
646     }
647
648     md->print_grpnms = NULL;
649
650     /* check whether we're going to write dh histograms */
651     md->dhc = NULL;
652     if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
653     {
654         /* Currently dh histograms are only written with dynamics */
655         if (EI_DYNAMICS(ir->eI))
656         {
657             snew(md->dhc, 1);
658
659             mde_delta_h_coll_init(md->dhc, ir);
660         }
661         md->fp_dhdl = NULL;
662         snew(md->dE, ir->fepvals->n_lambda);
663     }
664     else
665     {
666         md->fp_dhdl = fp_dhdl;
667         snew(md->dE, ir->fepvals->n_lambda);
668     }
669     if (ir->bSimTemp)
670     {
671         int i;
672         snew(md->temperatures, ir->fepvals->n_lambda);
673         for (i = 0; i < ir->fepvals->n_lambda; i++)
674         {
675             md->temperatures[i] = ir->simtempvals->temperatures[i];
676         }
677     }
678     return md;
679 }
680
681 /* print a lambda vector to a string
682    fep = the inputrec's FEP input data
683    i = the index of the lambda vector
684    get_native_lambda = whether to print the native lambda
685    get_names = whether to print the names rather than the values
686    str = the pre-allocated string buffer to print to. */
687 static void print_lambda_vector(t_lambda *fep, int i,
688                                 gmx_bool get_native_lambda, gmx_bool get_names,
689                                 char *str)
690 {
691     size_t nps = 0, np;
692     int    j, k = 0;
693     int    Nsep = 0;
694
695     for (j = 0; j < efptNR; j++)
696     {
697         if (fep->separate_dvdl[j])
698         {
699             Nsep++;
700         }
701     }
702     str[0] = 0; /* reset the string */
703     if (Nsep > 1)
704     {
705         str += sprintf(str, "("); /* set the opening parenthesis*/
706     }
707     for (j = 0; j < efptNR; j++)
708     {
709         if (fep->separate_dvdl[j])
710         {
711             double lam;
712             if (!get_names)
713             {
714                 if (get_native_lambda && fep->init_lambda >= 0)
715                 {
716                     str += sprintf(str, "%.4f", fep->init_lambda);
717                 }
718                 else
719                 {
720                     str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
721                 }
722             }
723             else
724             {
725                 str += sprintf(str, "%s", efpt_singular_names[j]);
726             }
727             /* print comma for the next item */
728             if (k < Nsep-1)
729             {
730                 str += sprintf(str, ", ");
731             }
732             k++;
733         }
734     }
735     if (Nsep > 1)
736     {
737         /* and add the closing parenthesis */
738         sprintf(str, ")");
739     }
740 }
741
742
743 extern FILE *open_dhdl(const char *filename, const t_inputrec *ir,
744                        const output_env_t oenv)
745 {
746     FILE       *fp;
747     const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
748     *lambdastate     = "\\lambda state", *remain = "remaining";
749     char        title[STRLEN], label_x[STRLEN], label_y[STRLEN];
750     int         i, np, nps, nsets, nsets_de, nsetsbegin;
751     int         n_lambda_terms = 0;
752     t_lambda   *fep            = ir->fepvals; /* for simplicity */
753     t_expanded *expand         = ir->expandedvals;
754     char      **setname;
755     char        buf[STRLEN], lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
756     int         bufplace = 0;
757
758     int         nsets_dhdl = 0;
759     int         s          = 0;
760     int         nsetsextend;
761     gmx_bool    write_pV = FALSE;
762
763     /* count the number of different lambda terms */
764     for (i = 0; i < efptNR; i++)
765     {
766         if (fep->separate_dvdl[i])
767         {
768             n_lambda_terms++;
769         }
770     }
771
772     if (fep->n_lambda == 0)
773     {
774         sprintf(title, "%s", dhdl);
775         sprintf(label_x, "Time (ps)");
776         sprintf(label_y, "%s (%s %s)",
777                 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
778     }
779     else
780     {
781         sprintf(title, "%s and %s", dhdl, deltag);
782         sprintf(label_x, "Time (ps)");
783         sprintf(label_y, "%s and %s (%s %s)",
784                 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
785     }
786     fp = gmx_fio_fopen(filename, "w+");
787     xvgr_header(fp, title, label_x, label_y, exvggtXNY, oenv);
788
789     if (!(ir->bSimTemp))
790     {
791         bufplace = sprintf(buf, "T = %g (K) ",
792                            ir->opts.ref_t[0]);
793     }
794     if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
795     {
796         if ( (fep->init_lambda >= 0)  && (n_lambda_terms == 1 ))
797         {
798             /* compatibility output */
799             sprintf(&(buf[bufplace]), "%s = %.4f", lambda, fep->init_lambda);
800         }
801         else
802         {
803             print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
804                                 lambda_vec_str);
805             print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
806                                 lambda_name_str);
807             sprintf(&(buf[bufplace]), "%s %d: %s = %s",
808                     lambdastate, fep->init_fep_state,
809                     lambda_name_str, lambda_vec_str);
810         }
811     }
812     xvgr_subtitle(fp, buf, oenv);
813
814
815     nsets_dhdl = 0;
816     if (fep->dhdl_derivatives == edhdlderivativesYES)
817     {
818         nsets_dhdl = n_lambda_terms;
819     }
820     /* count the number of delta_g states */
821     nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
822
823     nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
824
825     if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
826     {
827         nsets += 1;   /*add fep state for expanded ensemble */
828     }
829
830     if (fep->bPrintEnergy)
831     {
832         nsets += 1;  /* add energy to the dhdl as well */
833     }
834
835     nsetsextend = nsets;
836     if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
837     {
838         nsetsextend += 1; /* for PV term, other terms possible if required for
839                              the reduced potential (only needed with foreign
840                              lambda, and only output when init_lambda is not
841                              set in order to maintain compatibility of the
842                              dhdl.xvg file) */
843         write_pV     = TRUE;
844     }
845     snew(setname, nsetsextend);
846
847     if (expand->elmcmove > elmcmoveNO)
848     {
849         /* state for the fep_vals, if we have alchemical sampling */
850         sprintf(buf, "%s", "Thermodynamic state");
851         setname[s] = strdup(buf);
852         s         += 1;
853     }
854
855     if (fep->bPrintEnergy)
856     {
857         sprintf(buf, "%s (%s)", "Energy", unit_energy);
858         setname[s] = strdup(buf);
859         s         += 1;
860     }
861
862     if (fep->dhdl_derivatives == edhdlderivativesYES)
863     {
864         for (i = 0; i < efptNR; i++)
865         {
866             if (fep->separate_dvdl[i])
867             {
868
869                 if ( (fep->init_lambda >= 0)  && (n_lambda_terms == 1 ))
870                 {
871                     /* compatibility output */
872                     sprintf(buf, "%s %s %.4f", dhdl, lambda, fep->init_lambda);
873                 }
874                 else
875                 {
876                     double lam = fep->init_lambda;
877                     if (fep->init_lambda < 0)
878                     {
879                         lam = fep->all_lambda[i][fep->init_fep_state];
880                     }
881                     sprintf(buf, "%s %s = %.4f", dhdl, efpt_singular_names[i],
882                             lam);
883                 }
884                 setname[s] = strdup(buf);
885                 s         += 1;
886             }
887         }
888     }
889
890     if (fep->n_lambda > 0)
891     {
892         /* g_bar has to determine the lambda values used in this simulation
893          * from this xvg legend.
894          */
895
896         if (expand->elmcmove > elmcmoveNO)
897         {
898             nsetsbegin = 1;  /* for including the expanded ensemble */
899         }
900         else
901         {
902             nsetsbegin = 0;
903         }
904
905         if (fep->bPrintEnergy)
906         {
907             nsetsbegin += 1;
908         }
909         nsetsbegin += nsets_dhdl;
910
911         for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
912         {
913             print_lambda_vector(fep, i, FALSE, FALSE, lambda_vec_str);
914             if ( (fep->init_lambda >= 0)  && (n_lambda_terms == 1 ))
915             {
916                 /* for compatible dhdl.xvg files */
917                 nps = sprintf(buf, "%s %s %s", deltag, lambda, lambda_vec_str);
918             }
919             else
920             {
921                 nps = sprintf(buf, "%s %s to %s", deltag, lambda, lambda_vec_str);
922             }
923
924             if (ir->bSimTemp)
925             {
926                 /* print the temperature for this state if doing simulated annealing */
927                 sprintf(&buf[nps], "T = %g (%s)",
928                         ir->simtempvals->temperatures[s-(nsetsbegin)],
929                         unit_temp_K);
930             }
931             setname[s] = strdup(buf);
932             s++;
933         }
934         if (write_pV)
935         {
936             np                     = sprintf(buf, "pV (%s)", unit_energy);
937             setname[nsetsextend-1] = strdup(buf);  /* the first entry after
938                                                       nsets */
939         }
940
941         xvgr_legend(fp, nsetsextend, (const char **)setname, oenv);
942
943         for (s = 0; s < nsetsextend; s++)
944         {
945             sfree(setname[s]);
946         }
947         sfree(setname);
948     }
949
950     return fp;
951 }
952
953 static void copy_energy(t_mdebin *md, real e[], real ecpy[])
954 {
955     int i, j;
956
957     for (i = j = 0; (i < F_NRE); i++)
958     {
959         if (md->bEner[i])
960         {
961             ecpy[j++] = e[i];
962         }
963     }
964     if (j != md->f_nre)
965     {
966         gmx_incons("Number of energy terms wrong");
967     }
968 }
969
970 void upd_mdebin(t_mdebin       *md,
971                 gmx_bool        bDoDHDL,
972                 gmx_bool        bSum,
973                 double          time,
974                 real            tmass,
975                 gmx_enerdata_t *enerd,
976                 t_state        *state,
977                 t_lambda       *fep,
978                 t_expanded     *expand,
979                 matrix          box,
980                 tensor          svir,
981                 tensor          fvir,
982                 tensor          vir,
983                 tensor          pres,
984                 gmx_ekindata_t *ekind,
985                 rvec            mu_tot,
986                 gmx_constr_t    constr)
987 {
988     int    i, j, k, kk, m, n, gid;
989     real   crmsd[2], tmp6[6];
990     real   bs[NTRICLBOXS], vol, dens, pv, enthalpy;
991     real   eee[egNR];
992     real   ecopy[F_NRE];
993     double store_dhdl[efptNR];
994     real   store_energy = 0;
995     real   tmp;
996
997     /* Do NOT use the box in the state variable, but the separate box provided
998      * as an argument. This is because we sometimes need to write the box from
999      * the last timestep to match the trajectory frames.
1000      */
1001     copy_energy(md, enerd->term, ecopy);
1002     add_ebin(md->ebin, md->ie, md->f_nre, ecopy, bSum);
1003     if (md->nCrmsd)
1004     {
1005         crmsd[0] = constr_rmsd(constr, FALSE);
1006         if (md->nCrmsd > 1)
1007         {
1008             crmsd[1] = constr_rmsd(constr, TRUE);
1009         }
1010         add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
1011     }
1012     if (md->bDynBox)
1013     {
1014         int nboxs;
1015         if (md->bTricl)
1016         {
1017             bs[0] = box[XX][XX];
1018             bs[1] = box[YY][YY];
1019             bs[2] = box[ZZ][ZZ];
1020             bs[3] = box[YY][XX];
1021             bs[4] = box[ZZ][XX];
1022             bs[5] = box[ZZ][YY];
1023             nboxs = NTRICLBOXS;
1024         }
1025         else
1026         {
1027             bs[0] = box[XX][XX];
1028             bs[1] = box[YY][YY];
1029             bs[2] = box[ZZ][ZZ];
1030             nboxs = NBOXS;
1031         }
1032         vol  = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1033         dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1034         add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
1035         add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
1036         add_ebin(md->ebin, md->idens, 1, &dens, bSum);
1037
1038         if (md->bDiagPres)
1039         {
1040             /* This is pV (in kJ/mol).  The pressure is the reference pressure,
1041                not the instantaneous pressure */
1042             pv = vol*md->ref_p/PRESFAC;
1043
1044             add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
1045             enthalpy = pv + enerd->term[F_ETOT];
1046             add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
1047         }
1048     }
1049     if (md->bConstrVir)
1050     {
1051         add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
1052         add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
1053     }
1054     if (md->bVir)
1055     {
1056         add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
1057     }
1058     if (md->bPress)
1059     {
1060         add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
1061     }
1062     if (md->bSurft)
1063     {
1064         tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
1065         add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
1066     }
1067     if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
1068     {
1069         tmp6[0] = state->boxv[XX][XX];
1070         tmp6[1] = state->boxv[YY][YY];
1071         tmp6[2] = state->boxv[ZZ][ZZ];
1072         tmp6[3] = state->boxv[YY][XX];
1073         tmp6[4] = state->boxv[ZZ][XX];
1074         tmp6[5] = state->boxv[ZZ][YY];
1075         add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1076     }
1077     if (md->bMu)
1078     {
1079         add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1080     }
1081     if (ekind && ekind->cosacc.cos_accel != 0)
1082     {
1083         vol  = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1084         dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1085         add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1086         /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1087         tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1088                  *dens*sqr(box[ZZ][ZZ]*NANO/(2*M_PI)));
1089         add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1090     }
1091     if (md->nE > 1)
1092     {
1093         n = 0;
1094         for (i = 0; (i < md->nEg); i++)
1095         {
1096             for (j = i; (j < md->nEg); j++)
1097             {
1098                 gid = GID(i, j, md->nEg);
1099                 for (k = kk = 0; (k < egNR); k++)
1100                 {
1101                     if (md->bEInd[k])
1102                     {
1103                         eee[kk++] = enerd->grpp.ener[k][gid];
1104                     }
1105                 }
1106                 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1107                 n++;
1108             }
1109         }
1110     }
1111
1112     if (ekind)
1113     {
1114         for (i = 0; (i < md->nTC); i++)
1115         {
1116             md->tmp_r[i] = ekind->tcstat[i].T;
1117         }
1118         add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1119
1120         if (md->etc == etcNOSEHOOVER)
1121         {
1122             /* whether to print Nose-Hoover chains: */
1123             if (md->bPrintNHChains)
1124             {
1125                 if (md->bNHC_trotter)
1126                 {
1127                     for (i = 0; (i < md->nTC); i++)
1128                     {
1129                         for (j = 0; j < md->nNHC; j++)
1130                         {
1131                             k                = i*md->nNHC+j;
1132                             md->tmp_r[2*k]   = state->nosehoover_xi[k];
1133                             md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1134                         }
1135                     }
1136                     add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1137
1138                     if (md->bMTTK)
1139                     {
1140                         for (i = 0; (i < md->nTCP); i++)
1141                         {
1142                             for (j = 0; j < md->nNHC; j++)
1143                             {
1144                                 k                = i*md->nNHC+j;
1145                                 md->tmp_r[2*k]   = state->nhpres_xi[k];
1146                                 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1147                             }
1148                         }
1149                         add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1150                     }
1151                 }
1152                 else
1153                 {
1154                     for (i = 0; (i < md->nTC); i++)
1155                     {
1156                         md->tmp_r[2*i]   = state->nosehoover_xi[i];
1157                         md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1158                     }
1159                     add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1160                 }
1161             }
1162         }
1163         else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1164                  md->etc == etcVRESCALE)
1165         {
1166             for (i = 0; (i < md->nTC); i++)
1167             {
1168                 md->tmp_r[i] = ekind->tcstat[i].lambda;
1169             }
1170             add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1171         }
1172     }
1173
1174     if (ekind && md->nU > 1)
1175     {
1176         for (i = 0; (i < md->nU); i++)
1177         {
1178             copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1179         }
1180         add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1181     }
1182
1183     ebin_increase_count(md->ebin, bSum);
1184
1185     /* BAR + thermodynamic integration values */
1186     if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1187     {
1188         for (i = 0; i < enerd->n_lambda-1; i++)
1189         {
1190             /* zero for simulated tempering */
1191             md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1192             if (md->temperatures != NULL)
1193             {
1194                 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1195                 /* is this even useful to have at all? */
1196                 md->dE[i] += (md->temperatures[i]/
1197                               md->temperatures[state->fep_state]-1.0)*
1198                     enerd->term[F_EKIN];
1199             }
1200         }
1201
1202         if (md->fp_dhdl)
1203         {
1204             fprintf(md->fp_dhdl, "%.4f", time);
1205             /* the current free energy state */
1206
1207             /* print the current state if we are doing expanded ensemble */
1208             if (expand->elmcmove > elmcmoveNO)
1209             {
1210                 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1211             }
1212             /* total energy (for if the temperature changes */
1213             if (fep->bPrintEnergy)
1214             {
1215                 store_energy = enerd->term[F_ETOT];
1216                 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1217             }
1218
1219             if (fep->dhdl_derivatives == edhdlderivativesYES)
1220             {
1221                 for (i = 0; i < efptNR; i++)
1222                 {
1223                     if (fep->separate_dvdl[i])
1224                     {
1225                         /* assumes F_DVDL is first */
1226                         fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1227                     }
1228                 }
1229             }
1230             for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1231             {
1232                 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1233             }
1234             if ((md->epc != epcNO)  &&
1235                 (enerd->n_lambda > 0) &&
1236                 (fep->init_lambda < 0))
1237             {
1238                 fprintf(md->fp_dhdl, " %#.8g", pv);  /* PV term only needed when
1239                                                         there are alternate state
1240                                                         lambda and we're not in
1241                                                         compatibility mode */
1242             }
1243             fprintf(md->fp_dhdl, "\n");
1244             /* and the binary free energy output */
1245         }
1246         if (md->dhc && bDoDHDL)
1247         {
1248             int idhdl = 0;
1249             for (i = 0; i < efptNR; i++)
1250             {
1251                 if (fep->separate_dvdl[i])
1252                 {
1253                     /* assumes F_DVDL is first */
1254                     store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1255                     idhdl            += 1;
1256                 }
1257             }
1258             store_energy = enerd->term[F_ETOT];
1259             /* store_dh is dE */
1260             mde_delta_h_coll_add_dh(md->dhc,
1261                                     (double)state->fep_state,
1262                                     store_energy,
1263                                     pv,
1264                                     store_dhdl,
1265                                     md->dE + fep->lambda_start_n,
1266                                     time);
1267         }
1268     }
1269 }
1270
1271
1272 void upd_mdebin_step(t_mdebin *md)
1273 {
1274     ebin_increase_count(md->ebin, FALSE);
1275 }
1276
1277 static void npr(FILE *log, int n, char c)
1278 {
1279     for (; (n > 0); n--)
1280     {
1281         fprintf(log, "%c", c);
1282     }
1283 }
1284
1285 static void pprint(FILE *log, const char *s, t_mdebin *md)
1286 {
1287     char CHAR = '#';
1288     int  slen;
1289     char buf1[22], buf2[22];
1290
1291     slen = strlen(s);
1292     fprintf(log, "\t<======  ");
1293     npr(log, slen, CHAR);
1294     fprintf(log, "  ==>\n");
1295     fprintf(log, "\t<====  %s  ====>\n", s);
1296     fprintf(log, "\t<==  ");
1297     npr(log, slen, CHAR);
1298     fprintf(log, "  ======>\n\n");
1299
1300     fprintf(log, "\tStatistics over %s steps using %s frames\n",
1301             gmx_step_str(md->ebin->nsteps_sim, buf1),
1302             gmx_step_str(md->ebin->nsum_sim, buf2));
1303     fprintf(log, "\n");
1304 }
1305
1306 void print_ebin_header(FILE *log, gmx_int64_t steps, double time, real lambda)
1307 {
1308     char buf[22];
1309
1310     fprintf(log, "   %12s   %12s   %12s\n"
1311             "   %12s   %12.5f   %12.5f\n\n",
1312             "Step", "Time", "Lambda", gmx_step_str(steps, buf), time, lambda);
1313 }
1314
1315 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1316                 FILE *log,
1317                 gmx_int64_t step, double time,
1318                 int mode, gmx_bool bCompact,
1319                 t_mdebin *md, t_fcdata *fcd,
1320                 gmx_groups_t *groups, t_grpopts *opts)
1321 {
1322     /*static char **grpnms=NULL;*/
1323     char         buf[246];
1324     int          i, j, n, ni, nj, ndr, nor, b;
1325     int          ndisre = 0;
1326     real        *disre_rm3tav, *disre_rt;
1327
1328     /* these are for the old-style blocks (1 subblock, only reals), because
1329        there can be only one per ID for these */
1330     int          nr[enxNR];
1331     int          id[enxNR];
1332     real        *block[enxNR];
1333
1334     /* temporary arrays for the lambda values to write out */
1335     double      enxlambda_data[2];
1336
1337     t_enxframe  fr;
1338
1339     switch (mode)
1340     {
1341         case eprNORMAL:
1342             init_enxframe(&fr);
1343             fr.t            = time;
1344             fr.step         = step;
1345             fr.nsteps       = md->ebin->nsteps;
1346             fr.dt           = md->delta_t;
1347             fr.nsum         = md->ebin->nsum;
1348             fr.nre          = (bEne) ? md->ebin->nener : 0;
1349             fr.ener         = md->ebin->e;
1350             ndisre          = bDR ? fcd->disres.npair : 0;
1351             disre_rm3tav    = fcd->disres.rm3tav;
1352             disre_rt        = fcd->disres.rt;
1353             /* Optional additional old-style (real-only) blocks. */
1354             for (i = 0; i < enxNR; i++)
1355             {
1356                 nr[i] = 0;
1357             }
1358             if (fcd->orires.nr > 0 && bOR)
1359             {
1360                 diagonalize_orires_tensors(&(fcd->orires));
1361                 nr[enxOR]     = fcd->orires.nr;
1362                 block[enxOR]  = fcd->orires.otav;
1363                 id[enxOR]     = enxOR;
1364                 nr[enxORI]    = (fcd->orires.oinsl != fcd->orires.otav) ?
1365                     fcd->orires.nr : 0;
1366                 block[enxORI] = fcd->orires.oinsl;
1367                 id[enxORI]    = enxORI;
1368                 nr[enxORT]    = fcd->orires.nex*12;
1369                 block[enxORT] = fcd->orires.eig;
1370                 id[enxORT]    = enxORT;
1371             }
1372
1373             /* whether we are going to wrte anything out: */
1374             if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1375             {
1376
1377                 /* the old-style blocks go first */
1378                 fr.nblock = 0;
1379                 for (i = 0; i < enxNR; i++)
1380                 {
1381                     if (nr[i] > 0)
1382                     {
1383                         fr.nblock = i + 1;
1384                     }
1385                 }
1386                 add_blocks_enxframe(&fr, fr.nblock);
1387                 for (b = 0; b < fr.nblock; b++)
1388                 {
1389                     add_subblocks_enxblock(&(fr.block[b]), 1);
1390                     fr.block[b].id        = id[b];
1391                     fr.block[b].sub[0].nr = nr[b];
1392 #ifndef GMX_DOUBLE
1393                     fr.block[b].sub[0].type = xdr_datatype_float;
1394                     fr.block[b].sub[0].fval = block[b];
1395 #else
1396                     fr.block[b].sub[0].type = xdr_datatype_double;
1397                     fr.block[b].sub[0].dval = block[b];
1398 #endif
1399                 }
1400
1401                 /* check for disre block & fill it. */
1402                 if (ndisre > 0)
1403                 {
1404                     int db = fr.nblock;
1405                     fr.nblock += 1;
1406                     add_blocks_enxframe(&fr, fr.nblock);
1407
1408                     add_subblocks_enxblock(&(fr.block[db]), 2);
1409                     fr.block[db].id        = enxDISRE;
1410                     fr.block[db].sub[0].nr = ndisre;
1411                     fr.block[db].sub[1].nr = ndisre;
1412 #ifndef GMX_DOUBLE
1413                     fr.block[db].sub[0].type = xdr_datatype_float;
1414                     fr.block[db].sub[1].type = xdr_datatype_float;
1415                     fr.block[db].sub[0].fval = disre_rt;
1416                     fr.block[db].sub[1].fval = disre_rm3tav;
1417 #else
1418                     fr.block[db].sub[0].type = xdr_datatype_double;
1419                     fr.block[db].sub[1].type = xdr_datatype_double;
1420                     fr.block[db].sub[0].dval = disre_rt;
1421                     fr.block[db].sub[1].dval = disre_rm3tav;
1422 #endif
1423                 }
1424                 /* here we can put new-style blocks */
1425
1426                 /* Free energy perturbation blocks */
1427                 if (md->dhc)
1428                 {
1429                     mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1430                 }
1431
1432                 /* we can now free & reset the data in the blocks */
1433                 if (md->dhc)
1434                 {
1435                     mde_delta_h_coll_reset(md->dhc);
1436                 }
1437
1438                 /* do the actual I/O */
1439                 do_enx(fp_ene, &fr);
1440                 if (fr.nre)
1441                 {
1442                     /* We have stored the sums, so reset the sum history */
1443                     reset_ebin_sums(md->ebin);
1444                 }
1445             }
1446             free_enxframe(&fr);
1447             break;
1448         case eprAVER:
1449             if (log)
1450             {
1451                 pprint(log, "A V E R A G E S", md);
1452             }
1453             break;
1454         case eprRMS:
1455             if (log)
1456             {
1457                 pprint(log, "R M S - F L U C T U A T I O N S", md);
1458             }
1459             break;
1460         default:
1461             gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1462     }
1463
1464     if (log)
1465     {
1466         for (i = 0; i < opts->ngtc; i++)
1467         {
1468             if (opts->annealing[i] != eannNO)
1469             {
1470                 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1471                         *(groups->grpname[groups->grps[egcTC].nm_ind[i]]),
1472                         opts->ref_t[i]);
1473             }
1474         }
1475         if (mode == eprNORMAL && fcd->orires.nr > 0)
1476         {
1477             print_orires_log(log, &(fcd->orires));
1478         }
1479         fprintf(log, "   Energies (%s)\n", unit_energy);
1480         pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1481         fprintf(log, "\n");
1482
1483         if (!bCompact)
1484         {
1485             if (md->bDynBox)
1486             {
1487                 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1488                         mode, TRUE);
1489                 fprintf(log, "\n");
1490             }
1491             if (md->bConstrVir)
1492             {
1493                 fprintf(log, "   Constraint Virial (%s)\n", unit_energy);
1494                 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1495                 fprintf(log, "\n");
1496                 fprintf(log, "   Force Virial (%s)\n", unit_energy);
1497                 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1498                 fprintf(log, "\n");
1499             }
1500             if (md->bVir)
1501             {
1502                 fprintf(log, "   Total Virial (%s)\n", unit_energy);
1503                 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1504                 fprintf(log, "\n");
1505             }
1506             if (md->bPress)
1507             {
1508                 fprintf(log, "   Pressure (%s)\n", unit_pres_bar);
1509                 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1510                 fprintf(log, "\n");
1511             }
1512             if (md->bMu)
1513             {
1514                 fprintf(log, "   Total Dipole (%s)\n", unit_dipole_D);
1515                 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1516                 fprintf(log, "\n");
1517             }
1518
1519             if (md->nE > 1)
1520             {
1521                 if (md->print_grpnms == NULL)
1522                 {
1523                     snew(md->print_grpnms, md->nE);
1524                     n = 0;
1525                     for (i = 0; (i < md->nEg); i++)
1526                     {
1527                         ni = groups->grps[egcENER].nm_ind[i];
1528                         for (j = i; (j < md->nEg); j++)
1529                         {
1530                             nj = groups->grps[egcENER].nm_ind[j];
1531                             sprintf(buf, "%s-%s", *(groups->grpname[ni]),
1532                                     *(groups->grpname[nj]));
1533                             md->print_grpnms[n++] = strdup(buf);
1534                         }
1535                     }
1536                 }
1537                 sprintf(buf, "Epot (%s)", unit_energy);
1538                 fprintf(log, "%15s   ", buf);
1539                 for (i = 0; (i < egNR); i++)
1540                 {
1541                     if (md->bEInd[i])
1542                     {
1543                         fprintf(log, "%12s   ", egrp_nm[i]);
1544                     }
1545                 }
1546                 fprintf(log, "\n");
1547                 for (i = 0; (i < md->nE); i++)
1548                 {
1549                     fprintf(log, "%15s", md->print_grpnms[i]);
1550                     pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1551                             FALSE);
1552                 }
1553                 fprintf(log, "\n");
1554             }
1555             if (md->nTC > 1)
1556             {
1557                 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1558                 fprintf(log, "\n");
1559             }
1560             if (md->nU > 1)
1561             {
1562                 fprintf(log, "%15s   %12s   %12s   %12s\n",
1563                         "Group", "Ux", "Uy", "Uz");
1564                 for (i = 0; (i < md->nU); i++)
1565                 {
1566                     ni = groups->grps[egcACC].nm_ind[i];
1567                     fprintf(log, "%15s", *groups->grpname[ni]);
1568                     pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1569                 }
1570                 fprintf(log, "\n");
1571             }
1572         }
1573     }
1574
1575 }
1576
1577 void update_energyhistory(energyhistory_t * enerhist, t_mdebin * mdebin)
1578 {
1579     int i;
1580
1581     enerhist->nsteps     = mdebin->ebin->nsteps;
1582     enerhist->nsum       = mdebin->ebin->nsum;
1583     enerhist->nsteps_sim = mdebin->ebin->nsteps_sim;
1584     enerhist->nsum_sim   = mdebin->ebin->nsum_sim;
1585     enerhist->nener      = mdebin->ebin->nener;
1586
1587     if (mdebin->ebin->nsum > 0)
1588     {
1589         /* Check if we need to allocate first */
1590         if (enerhist->ener_ave == NULL)
1591         {
1592             snew(enerhist->ener_ave, enerhist->nener);
1593             snew(enerhist->ener_sum, enerhist->nener);
1594         }
1595
1596         for (i = 0; i < enerhist->nener; i++)
1597         {
1598             enerhist->ener_ave[i] = mdebin->ebin->e[i].eav;
1599             enerhist->ener_sum[i] = mdebin->ebin->e[i].esum;
1600         }
1601     }
1602
1603     if (mdebin->ebin->nsum_sim > 0)
1604     {
1605         /* Check if we need to allocate first */
1606         if (enerhist->ener_sum_sim == NULL)
1607         {
1608             snew(enerhist->ener_sum_sim, enerhist->nener);
1609         }
1610
1611         for (i = 0; i < enerhist->nener; i++)
1612         {
1613             enerhist->ener_sum_sim[i] = mdebin->ebin->e_sim[i].esum;
1614         }
1615     }
1616     if (mdebin->dhc)
1617     {
1618         mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1619     }
1620 }
1621
1622 void restore_energyhistory_from_state(t_mdebin        * mdebin,
1623                                       energyhistory_t * enerhist)
1624 {
1625     int i;
1626
1627     if ((enerhist->nsum > 0 || enerhist->nsum_sim > 0) &&
1628         mdebin->ebin->nener != enerhist->nener)
1629     {
1630         gmx_fatal(FARGS, "Mismatch between number of energies in run input (%d) and checkpoint file (%d).",
1631                   mdebin->ebin->nener, enerhist->nener);
1632     }
1633
1634     mdebin->ebin->nsteps     = enerhist->nsteps;
1635     mdebin->ebin->nsum       = enerhist->nsum;
1636     mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1637     mdebin->ebin->nsum_sim   = enerhist->nsum_sim;
1638
1639     for (i = 0; i < mdebin->ebin->nener; i++)
1640     {
1641         mdebin->ebin->e[i].eav  =
1642             (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1643         mdebin->ebin->e[i].esum =
1644             (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1645         mdebin->ebin->e_sim[i].esum =
1646             (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1647     }
1648     if (mdebin->dhc)
1649     {
1650         mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist);
1651     }
1652 }