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