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