2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/commandline/pargs.h"
46 #include "gromacs/commandline/viewit.h"
47 #include "gromacs/correlationfunctions/autocorr.h"
48 #include "gromacs/fileio/enxio.h"
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/fileio/tpxio.h"
51 #include "gromacs/fileio/trxio.h"
52 #include "gromacs/fileio/xvgr.h"
53 #include "gromacs/gmxana/gmx_ana.h"
54 #include "gromacs/gmxana/gstat.h"
55 #include "gromacs/math/functions.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdlib/mdebin.h"
59 #include "gromacs/mdtypes/inputrec.h"
60 #include "gromacs/mdtypes/md_enums.h"
61 #include "gromacs/topology/ifunc.h"
62 #include "gromacs/topology/mtop_lookup.h"
63 #include "gromacs/topology/mtop_util.h"
64 #include "gromacs/topology/topology.h"
65 #include "gromacs/trajectory/energyframe.h"
66 #include "gromacs/utility/arraysize.h"
67 #include "gromacs/utility/cstringutil.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/gmxassert.h"
70 #include "gromacs/utility/pleasecite.h"
71 #include "gromacs/utility/smalloc.h"
72 #include "gromacs/utility/strconvert.h"
74 static const int NOTSET
= -23451;
102 static void done_enerdata_t(int nset
, enerdata_t
*edat
)
107 for (int i
= 0; i
< nset
; i
++)
109 sfree(edat
->s
[i
].ener
);
110 sfree(edat
->s
[i
].es
);
115 static void chomp(char *buf
)
117 int len
= std::strlen(buf
);
119 while ((len
> 0) && (buf
[len
-1] == '\n'))
126 static int *select_by_name(int nre
, gmx_enxnm_t
*nm
, int *nset
)
129 int k
, kk
, j
, i
, nmatch
, nind
, nss
;
131 gmx_bool bEOF
, bVerbose
= TRUE
, bLong
= FALSE
;
132 char *ptr
, buf
[STRLEN
];
133 const char *fm4
= "%3d %-14s";
134 const char *fm2
= "%3d %-34s";
135 char **newnm
= nullptr;
137 if ((getenv("GMX_ENER_VERBOSE")) != nullptr)
142 fprintf(stderr
, "\n");
143 fprintf(stderr
, "Select the terms you want from the following list by\n");
144 fprintf(stderr
, "selecting either (part of) the name or the number or a combination.\n");
145 fprintf(stderr
, "End your selection with an empty line or a zero.\n");
146 fprintf(stderr
, "-------------------------------------------------------------------\n");
150 for (k
= 0; k
< nre
; k
++)
152 newnm
[k
] = gmx_strdup(nm
[k
].name
);
153 /* Insert dashes in all the names */
154 while ((ptr
= std::strchr(newnm
[k
], ' ')) != nullptr)
164 fprintf(stderr
, "\n");
167 for (kk
= k
; kk
< k
+4; kk
++)
169 if (kk
< nre
&& std::strlen(nm
[kk
].name
) > 14)
177 fprintf(stderr
, " ");
181 fprintf(stderr
, fm4
, k
+1, newnm
[k
]);
190 fprintf(stderr
, fm2
, k
+1, newnm
[k
]);
201 fprintf(stderr
, "\n\n");
207 while (!bEOF
&& (fgets2(buf
, STRLEN
-1, stdin
)))
209 /* Remove newlines */
215 /* Empty line means end of input */
216 bEOF
= (std::strlen(buf
) == 0);
224 /* First try to read an integer */
225 nss
= sscanf(ptr
, "%d", &nind
);
228 /* Zero means end of input */
233 else if ((1 <= nind
) && (nind
<= nre
))
239 fprintf(stderr
, "number %d is out of range\n", nind
);
244 /* Now try to read a string */
246 for (nind
= 0; nind
< nre
; nind
++)
248 if (gmx_strcasecmp(newnm
[nind
], ptr
) == 0)
256 i
= std::strlen(ptr
);
258 for (nind
= 0; nind
< nre
; nind
++)
260 if (gmx_strncasecmp(newnm
[nind
], ptr
, i
) == 0)
268 fprintf(stderr
, "String '%s' does not match anything\n", ptr
);
273 /* Look for the first space, and remove spaces from there */
274 if ((ptr
= std::strchr(ptr
, ' ')) != nullptr)
279 while (!bEOF
&& (ptr
&& (std::strlen(ptr
) > 0)));
284 for (i
= (*nset
) = 0; (i
< nre
); i
++)
296 gmx_fatal(FARGS
, "No energy terms selected");
299 for (i
= 0; (i
< nre
); i
++)
308 static void get_dhdl_parms(const char *topnm
, t_inputrec
*ir
)
314 /* all we need is the ir to be able to write the label */
315 read_tpx(topnm
, ir
, box
, &natoms
, nullptr, nullptr, &mtop
);
318 static void einstein_visco(const char *fn
, const char *fni
, int nsets
,
319 int nint
, real
**eneint
,
320 real V
, real T
, double dt
,
321 const gmx_output_env_t
*oenv
)
324 double av
[4], avold
[4];
330 for (i
= 0; i
<= nsets
; i
++)
334 fp0
= xvgropen(fni
, "Shear viscosity integral",
335 "Time (ps)", "(kg m\\S-1\\N s\\S-1\\N ps)", oenv
);
336 fp1
= xvgropen(fn
, "Shear viscosity using Einstein relation",
337 "Time (ps)", "(kg m\\S-1\\N s\\S-1\\N)", oenv
);
338 for (i
= 0; i
< nf4
; i
++)
340 for (m
= 0; m
<= nsets
; m
++)
344 for (j
= 0; j
< nint
- i
; j
++)
346 for (m
= 0; m
< nsets
; m
++)
348 di
= gmx::square(eneint
[m
][j
+i
] - eneint
[m
][j
]);
351 av
[nsets
] += di
/nsets
;
354 /* Convert to SI for the viscosity */
355 fac
= (V
*NANO
*NANO
*NANO
*PICO
*1e10
)/(2*BOLTZMANN
*T
)/(nint
- i
);
356 fprintf(fp0
, "%10g", i
*dt
);
357 for (m
= 0; (m
<= nsets
); m
++)
360 fprintf(fp0
, " %10g", av
[m
]);
363 fprintf(fp1
, "%10g", (i
+ 0.5)*dt
);
364 for (m
= 0; (m
<= nsets
); m
++)
366 fprintf(fp1
, " %10g", (av
[m
]-avold
[m
])/dt
);
389 static void clear_ee_sum(ee_sum_t
*ees
)
397 static void add_ee_sum(ee_sum_t
*ees
, double sum
, int np
)
403 static void add_ee_av(ee_sum_t
*ees
)
407 av
= ees
->sum
/ees
->np
;
414 static double calc_ee2(int nb
, ee_sum_t
*ees
)
416 return (ees
->sav2
/nb
- gmx::square(ees
->sav
/nb
))/(nb
- 1);
419 static void set_ee_av(ener_ee_t
*eee
)
423 char buf
[STEPSTRSIZE
];
424 fprintf(debug
, "Storing average for err.est.: %s steps\n",
425 gmx_step_str(eee
->nst
, buf
));
427 add_ee_av(&eee
->sum
);
429 if (eee
->b
== 1 || eee
->nst
< eee
->nst_min
)
431 eee
->nst_min
= eee
->nst
;
436 static void calc_averages(int nset
, enerdata_t
*edat
, int nbmin
, int nbmax
)
439 double sum
, sum2
, sump
, see2
;
440 int64_t np
, p
, bound_nb
;
444 double x
, sx
, sy
, sxx
, sxy
;
447 /* Check if we have exact statistics over all points */
448 for (i
= 0; i
< nset
; i
++)
451 ed
->bExactStat
= FALSE
;
454 /* All energy file sum entries 0 signals no exact sums.
455 * But if all energy values are 0, we still have exact sums.
458 for (f
= 0; f
< edat
->nframes
&& !ed
->bExactStat
; f
++)
460 if (ed
->ener
[i
] != 0)
464 ed
->bExactStat
= (ed
->es
[f
].sum
!= 0);
468 ed
->bExactStat
= TRUE
;
474 for (i
= 0; i
< nset
; i
++)
485 for (nb
= nbmin
; nb
<= nbmax
; nb
++)
488 clear_ee_sum(&eee
[nb
].sum
);
492 for (f
= 0; f
< edat
->nframes
; f
++)
498 /* Add the sum and the sum of variances to the totals. */
504 sum2
+= gmx::square(sum
/np
- (sum
+ es
->sum
)/(np
+ p
))
510 /* Add a single value to the sum and sum of squares. */
513 sum2
+= gmx::square(sump
);
516 /* sum has to be increased after sum2 */
520 /* For the linear regression use variance 1/p.
521 * Note that sump is the sum, not the average, so we don't need p*.
523 x
= edat
->step
[f
] - 0.5*(edat
->steps
[f
] - 1);
529 for (nb
= nbmin
; nb
<= nbmax
; nb
++)
531 /* Check if the current end step is closer to the desired
532 * block boundary than the next end step.
534 bound_nb
= (edat
->step
[0]-1)*nb
+ edat
->nsteps
*(eee
[nb
].b
+1);
535 if (eee
[nb
].nst
> 0 &&
536 bound_nb
- edat
->step
[f
-1]*nb
< edat
->step
[f
]*nb
- bound_nb
)
546 eee
[nb
].nst
+= edat
->step
[f
] - edat
->step
[f
-1];
550 add_ee_sum(&eee
[nb
].sum
, es
->sum
, edat
->points
[f
]);
554 add_ee_sum(&eee
[nb
].sum
, edat
->s
[i
].ener
[f
], 1);
556 bound_nb
= (edat
->step
[0]-1)*nb
+ edat
->nsteps
*(eee
[nb
].b
+1);
557 if (edat
->step
[f
]*nb
>= bound_nb
)
564 edat
->s
[i
].av
= sum
/np
;
567 edat
->s
[i
].rmsd
= std::sqrt(sum2
/np
);
571 edat
->s
[i
].rmsd
= std::sqrt(sum2
/np
- gmx::square(edat
->s
[i
].av
));
574 if (edat
->nframes
> 1)
576 edat
->s
[i
].slope
= (np
*sxy
- sx
*sy
)/(np
*sxx
- sx
*sx
);
580 edat
->s
[i
].slope
= 0;
585 for (nb
= nbmin
; nb
<= nbmax
; nb
++)
587 /* Check if we actually got nb blocks and if the smallest
588 * block is not shorter than 80% of the average.
592 char buf1
[STEPSTRSIZE
], buf2
[STEPSTRSIZE
];
593 fprintf(debug
, "Requested %d blocks, we have %d blocks, min %s nsteps %s\n",
595 gmx_step_str(eee
[nb
].nst_min
, buf1
),
596 gmx_step_str(edat
->nsteps
, buf2
));
598 if (eee
[nb
].b
== nb
&& 5*nb
*eee
[nb
].nst_min
>= 4*edat
->nsteps
)
600 see2
+= calc_ee2(nb
, &eee
[nb
].sum
);
606 edat
->s
[i
].ee
= std::sqrt(see2
/nee
);
616 static enerdata_t
*calc_sum(int nset
, enerdata_t
*edat
, int nbmin
, int nbmax
)
627 snew(s
->ener
, esum
->nframes
);
628 snew(s
->es
, esum
->nframes
);
630 s
->bExactStat
= TRUE
;
632 for (i
= 0; i
< nset
; i
++)
634 if (!edat
->s
[i
].bExactStat
)
636 s
->bExactStat
= FALSE
;
638 s
->slope
+= edat
->s
[i
].slope
;
641 for (f
= 0; f
< edat
->nframes
; f
++)
644 for (i
= 0; i
< nset
; i
++)
646 sum
+= edat
->s
[i
].ener
[f
];
650 for (i
= 0; i
< nset
; i
++)
652 sum
+= edat
->s
[i
].es
[f
].sum
;
658 calc_averages(1, esum
, nbmin
, nbmax
);
663 static void ee_pr(double ee
, int buflen
, char *buf
)
665 snprintf(buf
, buflen
, "%s", "--");
668 /* Round to two decimals by printing. */
670 snprintf(tmp
, sizeof(tmp
), "%.1e", ee
);
671 double rnd
= gmx::doubleFromString(tmp
);
672 snprintf(buf
, buflen
, "%g", rnd
);
676 static void remove_drift(int nset
, int nbmin
, int nbmax
, real dt
, enerdata_t
*edat
)
678 /* Remove the drift by performing a fit to y = ax+b.
679 Uses 5 iterations. */
683 edat
->npoints
= edat
->nframes
;
684 edat
->nsteps
= edat
->nframes
;
686 for (k
= 0; (k
< 5); k
++)
688 for (i
= 0; (i
< nset
); i
++)
690 delta
= edat
->s
[i
].slope
*dt
;
692 if (nullptr != debug
)
694 fprintf(debug
, "slope for set %d is %g\n", i
, edat
->s
[i
].slope
);
697 for (j
= 0; (j
< edat
->nframes
); j
++)
699 edat
->s
[i
].ener
[j
] -= j
*delta
;
700 edat
->s
[i
].es
[j
].sum
= 0;
701 edat
->s
[i
].es
[j
].sum2
= 0;
704 calc_averages(nset
, edat
, nbmin
, nbmax
);
708 static void calc_fluctuation_props(FILE *fp
,
709 gmx_bool bDriftCorr
, real dt
,
711 char **leg
, enerdata_t
*edat
,
712 int nbmin
, int nbmax
)
715 double vv
, v
, h
, varv
, hh
, varh
, tt
, cv
, cp
, alpha
, kappa
, dcp
, varet
;
718 eVol
, eEnth
, eTemp
, eEtot
, eNR
720 const char *my_ener
[] = { "Volume", "Enthalpy", "Temperature", "Total Energy" };
723 NANO3
= NANO
*NANO
*NANO
;
726 fprintf(fp
, "\nYou may want to use the -driftcorr flag in order to correct\n"
727 "for spurious drift in the graphs. Note that this is not\n"
728 "a substitute for proper equilibration and sampling!\n");
732 remove_drift(nset
, nbmin
, nbmax
, dt
, edat
);
734 for (i
= 0; (i
< eNR
); i
++)
736 for (ii
[i
] = 0; (ii
[i
] < nset
&&
737 (gmx_strcasecmp(leg
[ii
[i
]], my_ener
[i
]) != 0)); ii
[i
]++)
742 fprintf(fp,"Found %s data.\n",my_ener[i]);
744 /* Compute it all! */
745 alpha
= kappa
= cp
= dcp
= cv
= NOTSET
;
749 if (ii
[eTemp
] < nset
)
751 tt
= edat
->s
[ii
[eTemp
]].av
;
755 if ((ii
[eVol
] < nset
) && (ii
[eTemp
] < nset
))
757 vv
= edat
->s
[ii
[eVol
]].av
*NANO3
;
758 varv
= gmx::square(edat
->s
[ii
[eVol
]].rmsd
*NANO3
);
759 kappa
= (varv
/vv
)/(BOLTZMANN
*tt
);
763 if ((ii
[eEnth
] < nset
) && (ii
[eTemp
] < nset
))
765 hh
= KILO
*edat
->s
[ii
[eEnth
]].av
/AVOGADRO
;
766 varh
= gmx::square(KILO
*edat
->s
[ii
[eEnth
]].rmsd
/AVOGADRO
);
767 cp
= AVOGADRO
*((varh
/nmol
)/(BOLTZMANN
*tt
*tt
));
770 if ((ii
[eEtot
] < nset
) && (hh
== NOTSET
) && (tt
!= NOTSET
))
772 /* Only compute cv in constant volume runs, which we can test
773 by checking whether the enthalpy was computed.
775 varet
= gmx::square(edat
->s
[ii
[eEtot
]].rmsd
);
776 cv
= KILO
*((varet
/nmol
)/(BOLTZ
*tt
*tt
));
779 if ((ii
[eVol
] < nset
) && (ii
[eEnth
] < nset
) && (ii
[eTemp
] < nset
))
781 double v_sum
, h_sum
, vh_sum
, v_aver
, h_aver
, vh_aver
;
782 vh_sum
= v_sum
= h_sum
= 0;
783 for (j
= 0; (j
< edat
->nframes
); j
++)
785 v
= edat
->s
[ii
[eVol
]].ener
[j
]*NANO3
;
786 h
= KILO
*edat
->s
[ii
[eEnth
]].ener
[j
]/AVOGADRO
;
791 vh_aver
= vh_sum
/ edat
->nframes
;
792 v_aver
= v_sum
/ edat
->nframes
;
793 h_aver
= h_sum
/ edat
->nframes
;
794 alpha
= (vh_aver
-v_aver
*h_aver
)/(v_aver
*BOLTZMANN
*tt
*tt
);
795 dcp
= (v_aver
*AVOGADRO
/nmol
)*tt
*gmx::square(alpha
)/(kappa
);
802 fprintf(fp
, "\nWARNING: nmol = %d, this may not be what you want.\n",
805 fprintf(fp
, "\nTemperature dependent fluctuation properties at T = %g.\n", tt
);
806 fprintf(fp
, "\nHeat capacities obtained from fluctuations do *not* include\n");
807 fprintf(fp
, "quantum corrections. If you want to get a more accurate estimate\n");
808 fprintf(fp
, "please use the g_dos program.\n\n");
809 fprintf(fp
, "WARNING: Please verify that your simulations are converged and perform\n"
810 "a block-averaging error analysis (not implemented in g_energy yet)\n");
812 if (debug
!= nullptr)
816 fprintf(fp
, "varv = %10g (m^6)\n", varv
*AVOGADRO
/nmol
);
821 fprintf(fp
, "Volume = %10g m^3/mol\n",
826 fprintf(fp
, "Enthalpy = %10g kJ/mol\n",
827 hh
*AVOGADRO
/(KILO
*nmol
));
831 fprintf(fp
, "Coefficient of Thermal Expansion Alpha_P = %10g (1/K)\n",
836 fprintf(fp
, "Isothermal Compressibility Kappa = %10g (m^3/J)\n",
838 fprintf(fp
, "Adiabatic bulk modulus = %10g (J/m^3)\n",
843 fprintf(fp
, "Heat capacity at constant pressure Cp = %10g J/(mol K)\n",
848 fprintf(fp
, "Heat capacity at constant volume Cv = %10g J/(mol K)\n",
853 fprintf(fp
, "Cp-Cv = %10g J/(mol K)\n",
856 please_cite(fp
, "Allen1987a");
860 fprintf(fp
, "You should select the temperature in order to obtain fluctuation properties.\n");
864 static void analyse_ener(gmx_bool bCorr
, const char *corrfn
,
865 const char *eviscofn
, const char *eviscoifn
,
866 gmx_bool bFee
, gmx_bool bSum
, gmx_bool bFluct
,
867 gmx_bool bVisco
, const char *visfn
, int nmol
,
868 int64_t start_step
, double start_t
,
869 int64_t step
, double t
,
872 int nset
, const int set
[], const gmx_bool
*bIsEner
,
873 char **leg
, gmx_enxnm_t
*enm
,
874 real Vaver
, real ezero
,
875 int nbmin
, int nbmax
,
876 const gmx_output_env_t
*oenv
)
879 /* Check out the printed manual for equations! */
880 double Dt
, aver
, stddev
, errest
, delta_t
, totaldrift
;
881 enerdata_t
*esum
= nullptr;
882 real integral
, intBulk
, Temp
= 0, Pres
= 0;
883 real pr_aver
, pr_stddev
, pr_errest
;
884 double beta
= 0, expE
, expEtot
, *fee
= nullptr;
886 int nexact
, nnotexact
;
888 char buf
[256], eebuf
[100];
890 nsteps
= step
- start_step
+ 1;
893 fprintf(stdout
, "Not enough steps (%s) for statistics\n",
894 gmx_step_str(nsteps
, buf
));
898 /* Calculate the time difference */
899 delta_t
= t
- start_t
;
901 fprintf(stdout
, "\nStatistics over %s steps [ %.4f through %.4f ps ], %d data sets\n",
902 gmx_step_str(nsteps
, buf
), start_t
, t
, nset
);
904 calc_averages(nset
, edat
, nbmin
, nbmax
);
908 esum
= calc_sum(nset
, edat
, nbmin
, nbmax
);
911 if (!edat
->bHaveSums
)
920 for (i
= 0; (i
< nset
); i
++)
922 if (edat
->s
[i
].bExactStat
)
935 fprintf(stdout
, "All statistics are over %s points\n",
936 gmx_step_str(edat
->npoints
, buf
));
938 else if (nexact
== 0 || edat
->npoints
== edat
->nframes
)
940 fprintf(stdout
, "All statistics are over %d points (frames)\n",
945 fprintf(stdout
, "The term%s", nnotexact
== 1 ? "" : "s");
946 for (i
= 0; (i
< nset
); i
++)
948 if (!edat
->s
[i
].bExactStat
)
950 fprintf(stdout
, " '%s'", leg
[i
]);
953 fprintf(stdout
, " %s has statistics over %d points (frames)\n",
954 nnotexact
== 1 ? "is" : "are", edat
->nframes
);
955 fprintf(stdout
, "All other statistics are over %s points\n",
956 gmx_step_str(edat
->npoints
, buf
));
958 fprintf(stdout
, "\n");
960 fprintf(stdout
, "%-24s %10s %10s %10s %10s",
961 "Energy", "Average", "Err.Est.", "RMSD", "Tot-Drift");
964 fprintf(stdout
, " %10s\n", "-kT ln<e^(E/kT)>");
968 fprintf(stdout
, "\n");
970 fprintf(stdout
, "-------------------------------------------------------------------------------\n");
972 /* Initiate locals, only used with -sum */
976 beta
= 1.0/(BOLTZ
*reftemp
);
979 for (i
= 0; (i
< nset
); i
++)
981 aver
= edat
->s
[i
].av
;
982 stddev
= edat
->s
[i
].rmsd
;
983 errest
= edat
->s
[i
].ee
;
988 for (j
= 0; (j
< edat
->nframes
); j
++)
990 expE
+= std::exp(beta
*(edat
->s
[i
].ener
[j
] - aver
)/nmol
);
994 expEtot
+= expE
/edat
->nframes
;
997 fee
[i
] = std::log(expE
/edat
->nframes
)/beta
+ aver
/nmol
;
999 if (std::strstr(leg
[i
], "empera") != nullptr)
1003 else if (std::strstr(leg
[i
], "olum") != nullptr)
1007 else if (std::strstr(leg
[i
], "essure") != nullptr)
1013 pr_aver
= aver
/nmol
-ezero
;
1014 pr_stddev
= stddev
/nmol
;
1015 pr_errest
= errest
/nmol
;
1024 /* Multiply the slope in steps with the number of steps taken */
1025 totaldrift
= (edat
->nsteps
- 1)*edat
->s
[i
].slope
;
1031 ee_pr(pr_errest
, sizeof(eebuf
), eebuf
);
1032 fprintf(stdout
, "%-24s %10g %10s %10g %10g",
1033 leg
[i
], pr_aver
, eebuf
, pr_stddev
, totaldrift
);
1036 fprintf(stdout
, " %10g", fee
[i
]);
1039 fprintf(stdout
, " (%s)\n", enm
[set
[i
]].unit
);
1043 for (j
= 0; (j
< edat
->nframes
); j
++)
1045 edat
->s
[i
].ener
[j
] -= aver
;
1051 totaldrift
= (edat
->nsteps
- 1)*esum
->s
[0].slope
;
1052 ee_pr(esum
->s
[0].ee
/nmol
, sizeof(eebuf
), eebuf
);
1053 fprintf(stdout
, "%-24s %10g %10s %10s %10g (%s)",
1054 "Total", esum
->s
[0].av
/nmol
, eebuf
,
1055 "--", totaldrift
/nmol
, enm
[set
[0]].unit
);
1056 /* pr_aver,pr_stddev,a,totaldrift */
1059 fprintf(stdout
, " %10g %10g\n",
1060 std::log(expEtot
)/beta
+ esum
->s
[0].av
/nmol
, std::log(expEtot
)/beta
);
1064 fprintf(stdout
, "\n");
1068 /* Do correlation function */
1069 if (edat
->nframes
> 1)
1071 Dt
= delta_t
/(edat
->nframes
- 1);
1079 const char* leg
[] = { "Shear", "Bulk" };
1084 /* Assume pressure tensor is in Pxx Pxy Pxz Pyx Pyy Pyz Pzx Pzy Pzz */
1086 /* Symmetrise tensor! (and store in first three elements)
1087 * And subtract average pressure!
1090 for (i
= 0; i
< 12; i
++)
1092 snew(eneset
[i
], edat
->nframes
);
1094 for (i
= 0; (i
< edat
->nframes
); i
++)
1096 eneset
[0][i
] = 0.5*(edat
->s
[1].ener
[i
]+edat
->s
[3].ener
[i
]);
1097 eneset
[1][i
] = 0.5*(edat
->s
[2].ener
[i
]+edat
->s
[6].ener
[i
]);
1098 eneset
[2][i
] = 0.5*(edat
->s
[5].ener
[i
]+edat
->s
[7].ener
[i
]);
1099 for (j
= 3; j
<= 11; j
++)
1101 eneset
[j
][i
] = edat
->s
[j
].ener
[i
];
1103 eneset
[11][i
] -= Pres
;
1106 /* Determine integrals of the off-diagonal pressure elements */
1108 for (i
= 0; i
< 3; i
++)
1110 snew(eneint
[i
], edat
->nframes
+ 1);
1115 for (i
= 0; i
< edat
->nframes
; i
++)
1117 eneint
[0][i
+1] = eneint
[0][i
] + 0.5*(edat
->s
[1].es
[i
].sum
+ edat
->s
[3].es
[i
].sum
)*Dt
/edat
->points
[i
];
1118 eneint
[1][i
+1] = eneint
[1][i
] + 0.5*(edat
->s
[2].es
[i
].sum
+ edat
->s
[6].es
[i
].sum
)*Dt
/edat
->points
[i
];
1119 eneint
[2][i
+1] = eneint
[2][i
] + 0.5*(edat
->s
[5].es
[i
].sum
+ edat
->s
[7].es
[i
].sum
)*Dt
/edat
->points
[i
];
1122 einstein_visco(eviscofn
, eviscoifn
,
1123 3, edat
->nframes
+1, eneint
, Vaver
, Temp
, Dt
, oenv
);
1125 for (i
= 0; i
< 3; i
++)
1131 /*do_autocorr(corrfn,buf,nenergy,3,eneset,Dt,eacNormal,TRUE);*/
1132 /* Do it for shear viscosity */
1133 std::strcpy(buf
, "Shear Viscosity");
1134 low_do_autocorr(corrfn
, oenv
, buf
, edat
->nframes
, 3,
1135 (edat
->nframes
+1)/2, eneset
, Dt
,
1136 eacNormal
, 1, TRUE
, FALSE
, FALSE
, 0.0, 0.0, 0);
1138 /* Now for bulk viscosity */
1139 std::strcpy(buf
, "Bulk Viscosity");
1140 low_do_autocorr(corrfn
, oenv
, buf
, edat
->nframes
, 1,
1141 (edat
->nframes
+1)/2, &(eneset
[11]), Dt
,
1142 eacNormal
, 1, TRUE
, FALSE
, FALSE
, 0.0, 0.0, 0);
1144 factor
= (Vaver
*1e-26/(BOLTZMANN
*Temp
))*Dt
;
1145 fp
= xvgropen(visfn
, buf
, "Time (ps)", "\\8h\\4 (cp)", oenv
);
1146 xvgr_legend(fp
, asize(leg
), leg
, oenv
);
1148 /* Use trapezium rule for integration */
1151 nout
= get_acfnout();
1152 if ((nout
< 2) || (nout
>= edat
->nframes
/2))
1154 nout
= edat
->nframes
/2;
1156 for (i
= 1; (i
< nout
); i
++)
1158 integral
+= 0.5*(eneset
[0][i
-1] + eneset
[0][i
])*factor
;
1159 intBulk
+= 0.5*(eneset
[11][i
-1] + eneset
[11][i
])*factor
;
1160 fprintf(fp
, "%10g %10g %10g\n", (i
*Dt
), integral
, intBulk
);
1164 for (i
= 0; i
< 12; i
++)
1174 std::strcpy(buf
, "Autocorrelation of Energy Fluctuations");
1178 std::strcpy(buf
, "Energy Autocorrelation");
1181 do_autocorr(corrfn
, oenv
, buf
, edat
->nframes
,
1183 bSum
? &edat
->s
[nset
-1].ener
: eneset
,
1184 (delta_t
/edat
->nframes
), eacNormal
, FALSE
);
1190 static void print_time(FILE *fp
, double t
)
1192 fprintf(fp
, "%12.6f", t
);
1195 static void print1(FILE *fp
, gmx_bool bDp
, real e
)
1199 fprintf(fp
, " %16.12f", e
);
1203 fprintf(fp
, " %10.6f", e
);
1207 static void fec(const char *ene2fn
, const char *runavgfn
,
1208 real reftemp
, int nset
, const int set
[], char *leg
[],
1209 enerdata_t
*edat
, double time
[],
1210 const gmx_output_env_t
*oenv
)
1212 const char * ravgleg
[] = {
1213 "\\8D\\4E = E\\sB\\N-E\\sA\\N",
1214 "<e\\S-\\8D\\4E/kT\\N>\\s0..t\\N"
1218 int timecheck
, nenergy
, nenergy2
, maxenergy
;
1224 gmx_enxnm_t
*enm
= nullptr;
1228 /* read second energy file */
1231 enx
= open_enx(ene2fn
, "r");
1232 do_enxnms(enx
, &(fr
->nre
), &enm
);
1234 snew(eneset2
, nset
+1);
1240 /* This loop searches for the first frame (when -b option is given),
1241 * or when this has been found it reads just one energy frame
1245 bCont
= do_enx(enx
, fr
);
1249 timecheck
= check_times(fr
->t
);
1253 while (bCont
&& (timecheck
< 0));
1255 /* Store energies for analysis afterwards... */
1256 if ((timecheck
== 0) && bCont
)
1260 if (nenergy2
>= maxenergy
)
1263 for (i
= 0; i
<= nset
; i
++)
1265 srenew(eneset2
[i
], maxenergy
);
1268 GMX_RELEASE_ASSERT(time
!= nullptr, "trying to dereference NULL time pointer");
1270 if (fr
->t
!= time
[nenergy2
])
1272 fprintf(stderr
, "\nWARNING time mismatch %g!=%g at frame %s\n",
1273 fr
->t
, time
[nenergy2
], gmx_step_str(fr
->step
, buf
));
1275 for (i
= 0; i
< nset
; i
++)
1277 eneset2
[i
][nenergy2
] = fr
->ener
[set
[i
]].e
;
1283 while (bCont
&& (timecheck
== 0));
1286 if (edat
->nframes
!= nenergy2
)
1288 fprintf(stderr
, "\nWARNING file length mismatch %d!=%d\n",
1289 edat
->nframes
, nenergy2
);
1291 nenergy
= std::min(edat
->nframes
, nenergy2
);
1293 /* calculate fe difference dF = -kT ln < exp(-(E_B-E_A)/kT) >_A */
1297 fp
= xvgropen(runavgfn
, "Running average free energy difference",
1298 "Time (" unit_time
")", "\\8D\\4E (" unit_energy
")", oenv
);
1299 xvgr_legend(fp
, asize(ravgleg
), ravgleg
, oenv
);
1301 fprintf(stdout
, "\n%-24s %10s\n",
1302 "Energy", "dF = -kT ln < exp(-(EB-EA)/kT) >A");
1304 beta
= 1.0/(BOLTZ
*reftemp
);
1305 for (i
= 0; i
< nset
; i
++)
1307 if (gmx_strcasecmp(leg
[i
], enm
[set
[i
]].name
) != 0)
1309 fprintf(stderr
, "\nWARNING energy set name mismatch %s!=%s\n",
1310 leg
[i
], enm
[set
[i
]].name
);
1312 for (j
= 0; j
< nenergy
; j
++)
1314 dE
= eneset2
[i
][j
] - edat
->s
[i
].ener
[j
];
1315 sum
+= std::exp(-dE
*beta
);
1318 fprintf(fp
, "%10g %10g %10g\n",
1319 time
[j
], dE
, -BOLTZ
*reftemp
*std::log(sum
/(j
+1)) );
1322 aver
= -BOLTZ
*reftemp
*std::log(sum
/nenergy
);
1323 fprintf(stdout
, "%-24s %10g\n", leg
[i
], aver
);
1333 static void do_dhdl(t_enxframe
*fr
, const t_inputrec
*ir
, FILE **fp_dhdl
,
1334 const char *filename
, gmx_bool bDp
,
1335 int *blocks
, int *hists
, int *samples
, int *nlambdas
,
1336 const gmx_output_env_t
*oenv
)
1338 const char *dhdl
= "dH/d\\lambda", *deltag
= "\\DeltaH", *lambda
= "\\lambda";
1339 char title
[STRLEN
], label_x
[STRLEN
], label_y
[STRLEN
], legend
[STRLEN
];
1341 int nblock_hist
= 0, nblock_dh
= 0, nblock_dhcoll
= 0;
1344 double temp
= 0, start_time
= 0, delta_time
= 0, start_lambda
= 0;
1345 static int setnr
= 0;
1346 double *native_lambda_vec
= nullptr;
1347 const char **lambda_components
= nullptr;
1348 int n_lambda_vec
= 0;
1349 bool firstPass
= true;
1351 /* now count the blocks & handle the global dh data */
1352 for (i
= 0; i
< fr
->nblock
; i
++)
1354 if (fr
->block
[i
].id
== enxDHHIST
)
1358 else if (fr
->block
[i
].id
== enxDH
)
1362 else if (fr
->block
[i
].id
== enxDHCOLL
)
1365 if ( (fr
->block
[i
].nsub
< 1) ||
1366 (fr
->block
[i
].sub
[0].type
!= xdr_datatype_double
) ||
1367 (fr
->block
[i
].sub
[0].nr
< 5))
1369 gmx_fatal(FARGS
, "Unexpected block data");
1372 /* read the data from the DHCOLL block */
1373 temp
= fr
->block
[i
].sub
[0].dval
[0];
1374 start_time
= fr
->block
[i
].sub
[0].dval
[1];
1375 delta_time
= fr
->block
[i
].sub
[0].dval
[2];
1376 start_lambda
= fr
->block
[i
].sub
[0].dval
[3];
1377 if (fr
->block
[i
].nsub
> 1)
1381 n_lambda_vec
= fr
->block
[i
].sub
[1].ival
[1];
1382 snew(lambda_components
, n_lambda_vec
);
1383 snew(native_lambda_vec
, n_lambda_vec
);
1388 if (n_lambda_vec
!= fr
->block
[i
].sub
[1].ival
[1])
1391 "Unexpected change of basis set in lambda");
1394 for (j
= 0; j
< n_lambda_vec
; j
++)
1396 native_lambda_vec
[j
] = fr
->block
[i
].sub
[0].dval
[5+j
];
1397 lambda_components
[j
] =
1398 efpt_singular_names
[fr
->block
[i
].sub
[1].ival
[2+j
]];
1404 sfree(native_lambda_vec
);
1405 sfree(lambda_components
);
1406 if (nblock_hist
== 0 && nblock_dh
== 0)
1408 /* don't do anything */
1411 if (nblock_hist
> 0 && nblock_dh
> 0)
1413 gmx_fatal(FARGS
, "This energy file contains both histogram dhdl data and non-histogram dhdl data. Don't know what to do.");
1419 /* we have standard, non-histogram data --
1420 call open_dhdl to open the file */
1421 /* TODO this is an ugly hack that needs to be fixed: this will only
1422 work if the order of data is always the same and if we're
1423 only using the g_energy compiled with the mdrun that produced
1425 *fp_dhdl
= open_dhdl(filename
, ir
, oenv
);
1429 sprintf(title
, "N(%s)", deltag
);
1430 sprintf(label_x
, "%s (%s)", deltag
, unit_energy
);
1431 sprintf(label_y
, "Samples");
1432 *fp_dhdl
= xvgropen_type(filename
, title
, label_x
, label_y
, exvggtXNY
, oenv
);
1433 sprintf(buf
, "T = %g (K), %s = %g", temp
, lambda
, start_lambda
);
1434 xvgr_subtitle(*fp_dhdl
, buf
, oenv
);
1438 (*hists
) += nblock_hist
;
1439 (*blocks
) += nblock_dh
;
1440 (*nlambdas
) = nblock_hist
+nblock_dh
;
1442 /* write the data */
1443 if (nblock_hist
> 0)
1447 for (i
= 0; i
< fr
->nblock
; i
++)
1449 t_enxblock
*blk
= &(fr
->block
[i
]);
1450 if (blk
->id
== enxDHHIST
)
1452 double foreign_lambda
, dx
;
1454 int nhist
, derivative
;
1456 /* check the block types etc. */
1457 if ( (blk
->nsub
< 2) ||
1458 (blk
->sub
[0].type
!= xdr_datatype_double
) ||
1459 (blk
->sub
[1].type
!= xdr_datatype_int64
) ||
1460 (blk
->sub
[0].nr
< 2) ||
1461 (blk
->sub
[1].nr
< 2) )
1463 gmx_fatal(FARGS
, "Unexpected block data in file");
1465 foreign_lambda
= blk
->sub
[0].dval
[0];
1466 dx
= blk
->sub
[0].dval
[1];
1467 nhist
= blk
->sub
[1].lval
[0];
1468 derivative
= blk
->sub
[1].lval
[1];
1469 for (j
= 0; j
< nhist
; j
++)
1472 x0
= blk
->sub
[1].lval
[2+j
];
1476 sprintf(legend
, "N(%s(%s=%g) | %s=%g)",
1477 deltag
, lambda
, foreign_lambda
,
1478 lambda
, start_lambda
);
1482 sprintf(legend
, "N(%s | %s=%g)",
1483 dhdl
, lambda
, start_lambda
);
1487 xvgr_new_dataset(*fp_dhdl
, setnr
, 1, lg
, oenv
);
1489 for (k
= 0; k
< blk
->sub
[j
+2].nr
; k
++)
1494 hist
= blk
->sub
[j
+2].ival
[k
];
1497 fprintf(*fp_dhdl
, "%g %d\n%g %d\n", xmin
, hist
,
1501 /* multiple histogram data blocks in one histogram
1502 mean that the second one is the reverse of the first one:
1503 for dhdl derivatives, it's important to know both the
1504 maximum and minimum values */
1509 (*samples
) += static_cast<int>(sum
/nblock_hist
);
1516 for (i
= 0; i
< fr
->nblock
; i
++)
1518 t_enxblock
*blk
= &(fr
->block
[i
]);
1519 if (blk
->id
== enxDH
)
1523 len
= blk
->sub
[2].nr
;
1527 if (len
!= blk
->sub
[2].nr
)
1529 gmx_fatal(FARGS
, "Length inconsistency in dhdl data");
1536 for (i
= 0; i
< len
; i
++)
1538 double time
= start_time
+ delta_time
*i
;
1540 fprintf(*fp_dhdl
, "%.4f ", time
);
1542 for (j
= 0; j
< fr
->nblock
; j
++)
1544 t_enxblock
*blk
= &(fr
->block
[j
]);
1545 if (blk
->id
== enxDH
)
1548 if (blk
->sub
[2].type
== xdr_datatype_float
)
1550 value
= blk
->sub
[2].fval
[i
];
1554 value
= blk
->sub
[2].dval
[i
];
1556 /* we need to decide which data type it is based on the count*/
1558 if (j
== 1 && ir
->bExpanded
)
1560 fprintf(*fp_dhdl
, "%4d", static_cast<int>(value
)); /* if expanded ensembles and zero, this is a state value, it's an integer. We need a cleaner conditional than if j==1! */
1566 fprintf(*fp_dhdl
, " %#.12g", value
); /* print normal precision */
1570 fprintf(*fp_dhdl
, " %#.8g", value
); /* print normal precision */
1575 fprintf(*fp_dhdl
, "\n");
1581 int gmx_energy(int argc
, char *argv
[])
1583 const char *desc
[] = {
1584 "[THISMODULE] extracts energy components",
1585 "from an energy file. The user is prompted to interactively",
1586 "select the desired energy terms.[PAR]",
1588 "Average, RMSD, and drift are calculated with full precision from the",
1589 "simulation (see printed manual). Drift is calculated by performing",
1590 "a least-squares fit of the data to a straight line. The reported total drift",
1591 "is the difference of the fit at the first and last point.",
1592 "An error estimate of the average is given based on a block averages",
1593 "over 5 blocks using the full-precision averages. The error estimate",
1594 "can be performed over multiple block lengths with the options",
1595 "[TT]-nbmin[tt] and [TT]-nbmax[tt].",
1596 "[BB]Note[bb] that in most cases the energy files contains averages over all",
1597 "MD steps, or over many more points than the number of frames in",
1598 "energy file. This makes the [THISMODULE] statistics output more accurate",
1599 "than the [REF].xvg[ref] output. When exact averages are not present in the energy",
1600 "file, the statistics mentioned above are simply over the single, per-frame",
1601 "energy values.[PAR]",
1603 "The term fluctuation gives the RMSD around the least-squares fit.[PAR]",
1605 "Some fluctuation-dependent properties can be calculated provided",
1606 "the correct energy terms are selected, and that the command line option",
1607 "[TT]-fluct_props[tt] is given. The following properties",
1608 "will be computed:",
1610 "=============================== ===================",
1611 "Property Energy terms needed",
1612 "=============================== ===================",
1613 "Heat capacity C[SUB]p[sub] (NPT sims): Enthalpy, Temp",
1614 "Heat capacity C[SUB]v[sub] (NVT sims): Etot, Temp",
1615 "Thermal expansion coeff. (NPT): Enthalpy, Vol, Temp",
1616 "Isothermal compressibility: Vol, Temp",
1617 "Adiabatic bulk modulus: Vol, Temp",
1618 "=============================== ===================",
1620 "You always need to set the number of molecules [TT]-nmol[tt].",
1621 "The C[SUB]p[sub]/C[SUB]v[sub] computations do [BB]not[bb] include any corrections",
1622 "for quantum effects. Use the [gmx-dos] program if you need that (and you do).[PAR]",
1624 "Option [TT]-odh[tt] extracts and plots the free energy data",
1625 "(Hamiltoian differences and/or the Hamiltonian derivative dhdl)",
1626 "from the [TT]ener.edr[tt] file.[PAR]",
1628 "With [TT]-fee[tt] an estimate is calculated for the free-energy",
1629 "difference with an ideal gas state::",
1631 " [GRK]Delta[grk] A = A(N,V,T) - A[SUB]idealgas[sub](N,V,T) = kT [LN][CHEVRON][EXP]U[SUB]pot[sub]/kT[exp][chevron][ln]",
1632 " [GRK]Delta[grk] G = G(N,p,T) - G[SUB]idealgas[sub](N,p,T) = kT [LN][CHEVRON][EXP]U[SUB]pot[sub]/kT[exp][chevron][ln]",
1634 "where k is Boltzmann's constant, T is set by [TT]-fetemp[tt] and",
1635 "the average is over the ensemble (or time in a trajectory).",
1636 "Note that this is in principle",
1637 "only correct when averaging over the whole (Boltzmann) ensemble",
1638 "and using the potential energy. This also allows for an entropy",
1641 " [GRK]Delta[grk] S(N,V,T) = S(N,V,T) - S[SUB]idealgas[sub](N,V,T) = ([CHEVRON]U[SUB]pot[sub][chevron] - [GRK]Delta[grk] A)/T",
1642 " [GRK]Delta[grk] S(N,p,T) = S(N,p,T) - S[SUB]idealgas[sub](N,p,T) = ([CHEVRON]U[SUB]pot[sub][chevron] + pV - [GRK]Delta[grk] G)/T",
1645 "When a second energy file is specified ([TT]-f2[tt]), a free energy",
1646 "difference is calculated::",
1648 " dF = -kT [LN][CHEVRON][EXP]-(E[SUB]B[sub]-E[SUB]A[sub])/kT[exp][chevron][SUB]A[sub][ln] ,",
1650 "where E[SUB]A[sub] and E[SUB]B[sub] are the energies from the first and second energy",
1651 "files, and the average is over the ensemble A. The running average",
1652 "of the free energy difference is printed to a file specified by [TT]-ravg[tt].",
1653 "[BB]Note[bb] that the energies must both be calculated from the same trajectory."
1656 static gmx_bool bSum
= FALSE
, bFee
= FALSE
, bPrAll
= FALSE
, bFluct
= FALSE
, bDriftCorr
= FALSE
;
1657 static gmx_bool bDp
= FALSE
, bMutot
= FALSE
, bOrinst
= FALSE
, bOvec
= FALSE
, bFluctProps
= FALSE
;
1658 static int nmol
= 1, nbmin
= 5, nbmax
= 5;
1659 static real reftemp
= 300.0, ezero
= 0;
1661 { "-fee", FALSE
, etBOOL
, {&bFee
},
1662 "Do a free energy estimate" },
1663 { "-fetemp", FALSE
, etREAL
, {&reftemp
},
1664 "Reference temperature for free energy calculation" },
1665 { "-zero", FALSE
, etREAL
, {&ezero
},
1666 "Subtract a zero-point energy" },
1667 { "-sum", FALSE
, etBOOL
, {&bSum
},
1668 "Sum the energy terms selected rather than display them all" },
1669 { "-dp", FALSE
, etBOOL
, {&bDp
},
1670 "Print energies in high precision" },
1671 { "-nbmin", FALSE
, etINT
, {&nbmin
},
1672 "Minimum number of blocks for error estimate" },
1673 { "-nbmax", FALSE
, etINT
, {&nbmax
},
1674 "Maximum number of blocks for error estimate" },
1675 { "-mutot", FALSE
, etBOOL
, {&bMutot
},
1676 "Compute the total dipole moment from the components" },
1677 { "-aver", FALSE
, etBOOL
, {&bPrAll
},
1678 "Also print the exact average and rmsd stored in the energy frames (only when 1 term is requested)" },
1679 { "-nmol", FALSE
, etINT
, {&nmol
},
1680 "Number of molecules in your sample: the energies are divided by this number" },
1681 { "-fluct_props", FALSE
, etBOOL
, {&bFluctProps
},
1682 "Compute properties based on energy fluctuations, like heat capacity" },
1683 { "-driftcorr", FALSE
, etBOOL
, {&bDriftCorr
},
1684 "Useful only for calculations of fluctuation properties. The drift in the observables will be subtracted before computing the fluctuation properties."},
1685 { "-fluc", FALSE
, etBOOL
, {&bFluct
},
1686 "Calculate autocorrelation of energy fluctuations rather than energy itself" },
1687 { "-orinst", FALSE
, etBOOL
, {&bOrinst
},
1688 "Analyse instantaneous orientation data" },
1689 { "-ovec", FALSE
, etBOOL
, {&bOvec
},
1690 "Also plot the eigenvectors with [TT]-oten[tt]" }
1692 static const char *setnm
[] = {
1693 "Pres-XX", "Pres-XY", "Pres-XZ", "Pres-YX", "Pres-YY",
1694 "Pres-YZ", "Pres-ZX", "Pres-ZY", "Pres-ZZ", "Temperature",
1695 "Volume", "Pressure"
1698 FILE *out
= nullptr;
1699 FILE *fp_dhdl
= nullptr;
1703 gmx_enxnm_t
*enm
= nullptr;
1704 t_enxframe
*frame
, *fr
= nullptr;
1706 #define NEXT (1-cur)
1711 gmx_bool bFoundStart
, bCont
, bVisco
;
1713 double *time
= nullptr;
1715 int *set
= nullptr, i
, j
, nset
, sss
;
1716 gmx_bool
*bIsEner
= nullptr;
1717 char **leg
= nullptr;
1719 gmx_output_env_t
*oenv
;
1720 int dh_blocks
= 0, dh_hists
= 0, dh_samples
= 0, dh_lambdas
= 0;
1723 { efEDR
, "-f", nullptr, ffREAD
},
1724 { efEDR
, "-f2", nullptr, ffOPTRD
},
1725 { efTPR
, "-s", nullptr, ffOPTRD
},
1726 { efXVG
, "-o", "energy", ffWRITE
},
1727 { efXVG
, "-viol", "violaver", ffOPTWR
},
1728 { efXVG
, "-pairs", "pairs", ffOPTWR
},
1729 { efXVG
, "-corr", "enecorr", ffOPTWR
},
1730 { efXVG
, "-vis", "visco", ffOPTWR
},
1731 { efXVG
, "-evisco", "evisco", ffOPTWR
},
1732 { efXVG
, "-eviscoi", "eviscoi", ffOPTWR
},
1733 { efXVG
, "-ravg", "runavgdf", ffOPTWR
},
1734 { efXVG
, "-odh", "dhdl", ffOPTWR
}
1736 #define NFILE asize(fnm)
1741 ppa
= add_acf_pargs(&npargs
, pa
);
1742 if (!parse_common_args(&argc
, argv
,
1743 PCA_CAN_VIEW
| PCA_CAN_BEGIN
| PCA_CAN_END
,
1744 NFILE
, fnm
, npargs
, ppa
, asize(desc
), desc
, 0, nullptr, &oenv
))
1750 bDHDL
= opt2bSet("-odh", NFILE
, fnm
);
1755 fp
= open_enx(ftp2fn(efEDR
, NFILE
, fnm
), "r");
1756 do_enxnms(fp
, &nre
, &enm
);
1760 bVisco
= opt2bSet("-vis", NFILE
, fnm
);
1762 t_inputrec irInstance
;
1763 t_inputrec
*ir
= &irInstance
;
1769 nset
= asize(setnm
);
1771 /* This is nasty code... To extract Pres tensor, Volume and Temperature */
1772 for (j
= 0; j
< nset
; j
++)
1774 for (i
= 0; i
< nre
; i
++)
1776 if (std::strstr(enm
[i
].name
, setnm
[j
]))
1784 if (gmx_strcasecmp(setnm
[j
], "Volume") == 0)
1786 printf("Enter the box volume (" unit_volume
"): ");
1787 if (1 != scanf("%lf", &dbl
))
1789 gmx_fatal(FARGS
, "Error reading user input");
1795 gmx_fatal(FARGS
, "Could not find term %s for viscosity calculation",
1803 set
= select_by_name(nre
, enm
, &nset
);
1805 /* Print all the different units once */
1806 sprintf(buf
, "(%s)", enm
[set
[0]].unit
);
1807 for (i
= 1; i
< nset
; i
++)
1809 for (j
= 0; j
< i
; j
++)
1811 if (std::strcmp(enm
[set
[i
]].unit
, enm
[set
[j
]].unit
) == 0)
1818 std::strcat(buf
, ", (");
1819 std::strcat(buf
, enm
[set
[i
]].unit
);
1820 std::strcat(buf
, ")");
1823 out
= xvgropen(opt2fn("-o", NFILE
, fnm
), "GROMACS Energies", "Time (ps)", buf
,
1827 for (i
= 0; (i
< nset
); i
++)
1829 leg
[i
] = enm
[set
[i
]].name
;
1833 leg
[nset
] = gmx_strdup("Sum");
1834 xvgr_legend(out
, nset
+1, (const char**)leg
, oenv
);
1838 xvgr_legend(out
, nset
, (const char**)leg
, oenv
);
1841 snew(bIsEner
, nset
);
1842 for (i
= 0; (i
< nset
); i
++)
1845 for (j
= 0; (j
<= F_ETOT
); j
++)
1847 bIsEner
[i
] = bIsEner
[i
] ||
1848 (gmx_strcasecmp(interaction_function
[j
].longname
, leg
[i
]) == 0);
1851 if (bPrAll
&& nset
> 1)
1853 gmx_fatal(FARGS
, "Printing averages can only be done when a single set is selected");
1859 get_dhdl_parms(ftp2fn(efTPR
, NFILE
, fnm
), ir
);
1862 /* Initiate energies and set them to zero */
1866 edat
.step
= nullptr;
1867 edat
.steps
= nullptr;
1868 edat
.points
= nullptr;
1869 edat
.bHaveSums
= TRUE
;
1872 /* Initiate counters */
1873 bFoundStart
= FALSE
;
1878 /* This loop searches for the first frame (when -b option is given),
1879 * or when this has been found it reads just one energy frame
1883 bCont
= do_enx(fp
, &(frame
[NEXT
]));
1886 timecheck
= check_times(frame
[NEXT
].t
);
1889 while (bCont
&& (timecheck
< 0));
1891 if ((timecheck
== 0) && bCont
)
1893 /* We read a valid frame, so we can use it */
1894 fr
= &(frame
[NEXT
]);
1898 /* The frame contains energies, so update cur */
1901 if (edat
.nframes
% 1000 == 0)
1903 srenew(edat
.step
, edat
.nframes
+1000);
1904 std::memset(&(edat
.step
[edat
.nframes
]), 0, 1000*sizeof(edat
.step
[0]));
1905 srenew(edat
.steps
, edat
.nframes
+1000);
1906 std::memset(&(edat
.steps
[edat
.nframes
]), 0, 1000*sizeof(edat
.steps
[0]));
1907 srenew(edat
.points
, edat
.nframes
+1000);
1908 std::memset(&(edat
.points
[edat
.nframes
]), 0, 1000*sizeof(edat
.points
[0]));
1910 for (i
= 0; i
< nset
; i
++)
1912 srenew(edat
.s
[i
].ener
, edat
.nframes
+1000);
1913 std::memset(&(edat
.s
[i
].ener
[edat
.nframes
]), 0,
1914 1000*sizeof(edat
.s
[i
].ener
[0]));
1915 srenew(edat
.s
[i
].es
, edat
.nframes
+1000);
1916 std::memset(&(edat
.s
[i
].es
[edat
.nframes
]), 0,
1917 1000*sizeof(edat
.s
[i
].es
[0]));
1922 edat
.step
[nfr
] = fr
->step
;
1927 /* Initiate the previous step data */
1928 start_step
= fr
->step
;
1930 /* Initiate the energy sums */
1931 edat
.steps
[nfr
] = 1;
1932 edat
.points
[nfr
] = 1;
1933 for (i
= 0; i
< nset
; i
++)
1936 edat
.s
[i
].es
[nfr
].sum
= fr
->ener
[sss
].e
;
1937 edat
.s
[i
].es
[nfr
].sum2
= 0;
1944 edat
.steps
[nfr
] = fr
->nsteps
;
1948 /* mdrun only calculated the energy at energy output
1949 * steps. We don't need to check step intervals.
1951 edat
.points
[nfr
] = 1;
1952 for (i
= 0; i
< nset
; i
++)
1955 edat
.s
[i
].es
[nfr
].sum
= fr
->ener
[sss
].e
;
1956 edat
.s
[i
].es
[nfr
].sum2
= 0;
1959 edat
.bHaveSums
= FALSE
;
1961 else if (fr
->step
- start_step
+ 1 == edat
.nsteps
+ fr
->nsteps
)
1963 /* We have statistics to the previous frame */
1964 edat
.points
[nfr
] = fr
->nsum
;
1965 for (i
= 0; i
< nset
; i
++)
1968 edat
.s
[i
].es
[nfr
].sum
= fr
->ener
[sss
].esum
;
1969 edat
.s
[i
].es
[nfr
].sum2
= fr
->ener
[sss
].eav
;
1971 edat
.npoints
+= fr
->nsum
;
1975 /* The interval does not match fr->nsteps:
1976 * can not do exact averages.
1978 edat
.bHaveSums
= FALSE
;
1981 edat
.nsteps
= fr
->step
- start_step
+ 1;
1983 for (i
= 0; i
< nset
; i
++)
1985 edat
.s
[i
].ener
[nfr
] = fr
->ener
[set
[i
]].e
;
1989 * Store energies for analysis afterwards...
1991 if (!bDHDL
&& (fr
->nre
> 0))
1993 if (edat
.nframes
% 1000 == 0)
1995 srenew(time
, edat
.nframes
+1000);
1997 time
[edat
.nframes
] = fr
->t
;
2002 do_dhdl(fr
, ir
, &fp_dhdl
, opt2fn("-odh", NFILE
, fnm
), bDp
, &dh_blocks
, &dh_hists
, &dh_samples
, &dh_lambdas
, oenv
);
2005 /*******************************************
2007 *******************************************/
2014 /* We skip frames with single points (usually only the first frame),
2015 * since they would result in an average plot with outliers.
2019 print_time(out
, fr
->t
);
2020 print1(out
, bDp
, fr
->ener
[set
[0]].e
);
2021 print1(out
, bDp
, fr
->ener
[set
[0]].esum
/fr
->nsum
);
2022 print1(out
, bDp
, std::sqrt(fr
->ener
[set
[0]].eav
/fr
->nsum
));
2028 print_time(out
, fr
->t
);
2032 for (i
= 0; i
< nset
; i
++)
2034 sum
+= fr
->ener
[set
[i
]].e
;
2036 print1(out
, bDp
, sum
/nmol
-ezero
);
2040 for (i
= 0; (i
< nset
); i
++)
2044 print1(out
, bDp
, (fr
->ener
[set
[i
]].e
)/nmol
-ezero
);
2048 print1(out
, bDp
, fr
->ener
[set
[i
]].e
);
2058 while (bCont
&& (timecheck
== 0));
2060 fprintf(stderr
, "\n");
2071 gmx_fio_fclose(fp_dhdl
);
2072 printf("\n\nWrote %d lambda values with %d samples as ",
2073 dh_lambdas
, dh_samples
);
2076 printf("%d dH histograms ", dh_hists
);
2080 printf("%d dH data blocks ", dh_blocks
);
2082 printf("to %s\n", opt2fn("-odh", NFILE
, fnm
));
2087 gmx_fatal(FARGS
, "No dH data in %s\n", opt2fn("-f", NFILE
, fnm
));
2093 double dt
= (frame
[cur
].t
-start_t
)/(edat
.nframes
-1);
2094 analyse_ener(opt2bSet("-corr", NFILE
, fnm
), opt2fn("-corr", NFILE
, fnm
),
2095 opt2fn("-evisco", NFILE
, fnm
), opt2fn("-eviscoi", NFILE
, fnm
),
2097 bVisco
, opt2fn("-vis", NFILE
, fnm
),
2099 start_step
, start_t
, frame
[cur
].step
, frame
[cur
].t
,
2101 nset
, set
, bIsEner
, leg
, enm
, Vaver
, ezero
, nbmin
, nbmax
,
2105 calc_fluctuation_props(stdout
, bDriftCorr
, dt
, nset
, nmol
, leg
, &edat
,
2109 if (opt2bSet("-f2", NFILE
, fnm
))
2111 fec(opt2fn("-f2", NFILE
, fnm
), opt2fn("-ravg", NFILE
, fnm
),
2112 reftemp
, nset
, set
, leg
, &edat
, time
, oenv
);
2115 done_enerdata_t(nset
, &edat
);
2117 free_enxframe(&frame
[0]);
2118 free_enxframe(&frame
[1]);
2120 free_enxnms(nre
, enm
);
2126 const char *nxy
= "-nxy";
2128 do_view(oenv
, opt2fn("-o", NFILE
, fnm
), nxy
);
2129 do_view(oenv
, opt2fn_null("-ravg", NFILE
, fnm
), nxy
);
2130 do_view(oenv
, opt2fn_null("-odh", NFILE
, fnm
), nxy
);
2132 output_env_done(oenv
);