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