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