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/awh/awh.h"
46 #include "gromacs/fileio/enxio.h"
47 #include "gromacs/fileio/gmxfio.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/listed-forces/disre.h"
51 #include "gromacs/listed-forces/orires.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdlib/constr.h"
56 #include "gromacs/mdlib/mdebin_bar.h"
57 #include "gromacs/mdlib/mdrun.h"
58 #include "gromacs/mdtypes/energyhistory.h"
59 #include "gromacs/mdtypes/fcdata.h"
60 #include "gromacs/mdtypes/group.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/mdtypes/state.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/pulling/pull.h"
66 #include "gromacs/topology/mtop_util.h"
67 #include "gromacs/trajectory/energyframe.h"
68 #include "gromacs/utility/arraysize.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/smalloc.h"
72 static const char *conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
74 static const char *boxs_nm[] = { "Box-X", "Box-Y", "Box-Z" };
76 static const char *tricl_boxs_nm[] = {
77 "Box-XX", "Box-YY", "Box-ZZ",
78 "Box-YX", "Box-ZX", "Box-ZY"
81 static const char *vol_nm[] = { "Volume" };
83 static const char *dens_nm[] = {"Density" };
85 static const char *pv_nm[] = {"pV" };
87 static const char *enthalpy_nm[] = {"Enthalpy" };
89 static const char *boxvel_nm[] = {
90 "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
91 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
94 #define NBOXS asize(boxs_nm)
95 #define NTRICLBOXS asize(tricl_boxs_nm)
97 t_mdebin *init_mdebin(ener_file_t fp_ene,
98 const gmx_mtop_t *mtop,
102 const char *ener_nm[F_NRE];
103 static const char *vir_nm[] = {
104 "Vir-XX", "Vir-XY", "Vir-XZ",
105 "Vir-YX", "Vir-YY", "Vir-YZ",
106 "Vir-ZX", "Vir-ZY", "Vir-ZZ"
108 static const char *sv_nm[] = {
109 "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
110 "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
111 "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ"
113 static const char *fv_nm[] = {
114 "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
115 "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
116 "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ"
118 static const char *pres_nm[] = {
119 "Pres-XX", "Pres-XY", "Pres-XZ",
120 "Pres-YX", "Pres-YY", "Pres-YZ",
121 "Pres-ZX", "Pres-ZY", "Pres-ZZ"
123 static const char *surft_nm[] = {
126 static const char *mu_nm[] = {
127 "Mu-X", "Mu-Y", "Mu-Z"
129 static const char *vcos_nm[] = {
132 static const char *visc_nm[] = {
135 static const char *baro_nm[] = {
139 const gmx_groups_t *groups;
144 int i, j, ni, nj, n, k, kk, ncon, nset;
149 if (EI_DYNAMICS(ir->eI))
151 md->delta_t = ir->delta_t;
158 groups = &mtop->groups;
160 bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
161 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
162 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
164 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
165 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
166 md->bConstr = (ncon > 0 || nset > 0);
167 md->bConstrVir = FALSE;
170 if (ncon > 0 && ir->eConstrAlg == econtLINCS)
174 md->bConstrVir = (getenv("GMX_CONSTRAINTVIR") != nullptr);
181 /* Energy monitoring */
182 for (i = 0; i < egNR; i++)
184 md->bEInd[i] = FALSE;
187 for (i = 0; i < F_NRE; i++)
189 md->bEner[i] = FALSE;
192 md->bEner[i] = !bBHAM;
194 else if (i == F_BHAM)
196 md->bEner[i] = bBHAM;
200 md->bEner[i] = ir->bQMMM;
202 else if (i == F_RF_EXCL)
204 md->bEner[i] = (EEL_RF(ir->coulombtype) && ir->cutoff_scheme == ecutsGROUP);
206 else if (i == F_COUL_RECIP)
208 md->bEner[i] = EEL_FULL(ir->coulombtype);
210 else if (i == F_LJ_RECIP)
212 md->bEner[i] = EVDW_PME(ir->vdwtype);
214 else if (i == F_LJ14)
218 else if (i == F_COUL14)
222 else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
224 md->bEner[i] = FALSE;
226 else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
227 (i == F_DVDL_VDW && ir->fepvals->separate_dvdl[efptVDW]) ||
228 (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
229 (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
230 (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
231 (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
233 md->bEner[i] = (ir->efep != efepNO);
235 else if ((interaction_function[i].flags & IF_VSITE) ||
236 (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
238 md->bEner[i] = FALSE;
240 else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES) || (i == F_EQM))
244 else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
246 md->bEner[i] = EI_DYNAMICS(ir->eI);
248 else if (i == F_DISPCORR || i == F_PDISPCORR)
250 md->bEner[i] = (ir->eDispCorr != edispcNO);
252 else if (i == F_DISRESVIOL)
254 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
256 else if (i == F_ORIRESDEV)
258 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
260 else if (i == F_CONNBONDS)
262 md->bEner[i] = FALSE;
264 else if (i == F_COM_PULL)
266 md->bEner[i] = (ir->bPull && pull_have_potential(ir->pull_work));
268 else if (i == F_ECONSERVED)
270 md->bEner[i] = (integratorHasConservedEnergyQuantity(ir));
274 md->bEner[i] = (gmx_mtop_ftype_count(mtop, i) > 0);
279 for (i = 0; i < F_NRE; i++)
283 ener_nm[md->f_nre] = interaction_function[i].longname;
289 md->bDiagPres = !TRICLINIC(ir->ref_p);
290 md->ref_p = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
291 md->bTricl = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
292 md->bDynBox = inputrecDynamicBox(ir);
294 md->bNHC_trotter = inputrecNvtTrotter(ir);
295 md->bPrintNHChains = ir->bPrintNHChains;
296 md->bMTTK = (inputrecNptTrotter(ir) || inputrecNphTrotter(ir));
297 md->bMu = inputrecNeedMutot(ir);
299 md->ebin = mk_ebin();
300 /* Pass NULL for unit to let get_ebin_space determine the units
301 * for interaction_function[i].longname
303 md->ie = get_ebin_space(md->ebin, md->f_nre, ener_nm, nullptr);
306 /* This should be called directly after the call for md->ie,
307 * such that md->iconrmsd follows directly in the list.
309 md->iconrmsd = get_ebin_space(md->ebin, md->nCrmsd, conrmsd_nm, "");
313 md->ib = get_ebin_space(md->ebin,
314 md->bTricl ? NTRICLBOXS : NBOXS,
315 md->bTricl ? tricl_boxs_nm : boxs_nm,
317 md->ivol = get_ebin_space(md->ebin, 1, vol_nm, unit_volume);
318 md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
321 md->ipv = get_ebin_space(md->ebin, 1, pv_nm, unit_energy);
322 md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm, unit_energy);
327 md->isvir = get_ebin_space(md->ebin, asize(sv_nm), sv_nm, unit_energy);
328 md->ifvir = get_ebin_space(md->ebin, asize(fv_nm), fv_nm, unit_energy);
330 md->ivir = get_ebin_space(md->ebin, asize(vir_nm), vir_nm, unit_energy);
331 md->ipres = get_ebin_space(md->ebin, asize(pres_nm), pres_nm, unit_pres_bar);
332 md->isurft = get_ebin_space(md->ebin, asize(surft_nm), surft_nm,
334 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
336 md->ipc = get_ebin_space(md->ebin, md->bTricl ? 6 : 3,
337 boxvel_nm, unit_vel);
341 md->imu = get_ebin_space(md->ebin, asize(mu_nm), mu_nm, unit_dipole_D);
343 if (ir->cos_accel != 0)
345 md->ivcos = get_ebin_space(md->ebin, asize(vcos_nm), vcos_nm, unit_vel);
346 md->ivisc = get_ebin_space(md->ebin, asize(visc_nm), visc_nm,
350 /* Energy monitoring */
351 for (i = 0; i < egNR; i++)
353 md->bEInd[i] = FALSE;
355 md->bEInd[egCOULSR] = TRUE;
356 md->bEInd[egLJSR ] = TRUE;
360 md->bEInd[egLJSR] = FALSE;
361 md->bEInd[egBHAMSR] = TRUE;
365 md->bEInd[egLJ14] = TRUE;
366 md->bEInd[egCOUL14] = TRUE;
369 for (i = 0; (i < egNR); i++)
377 n = groups->grps[egcENER].nr;
379 md->nE = (n*(n+1))/2;
381 snew(md->igrp, md->nE);
386 for (k = 0; (k < md->nEc); k++)
388 snew(gnm[k], STRLEN);
390 for (i = 0; (i < groups->grps[egcENER].nr); i++)
392 ni = groups->grps[egcENER].nm_ind[i];
393 for (j = i; (j < groups->grps[egcENER].nr); j++)
395 nj = groups->grps[egcENER].nm_ind[j];
396 for (k = kk = 0; (k < egNR); k++)
400 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
401 *(groups->grpname[ni]), *(groups->grpname[nj]));
405 md->igrp[n] = get_ebin_space(md->ebin, md->nEc,
406 (const char **)gnm, unit_energy);
410 for (k = 0; (k < md->nEc); k++)
418 gmx_incons("Number of energy terms wrong");
422 md->nTC = groups->grps[egcTC].nr;
423 md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
426 md->nTCP = 1; /* assume only one possible coupling system for barostat
433 if (md->etc == etcNOSEHOOVER)
435 if (md->bNHC_trotter)
437 md->mde_n = 2*md->nNHC*md->nTC;
441 md->mde_n = 2*md->nTC;
443 if (md->epc == epcMTTK)
445 md->mdeb_n = 2*md->nNHC*md->nTCP;
454 snew(md->tmp_r, md->mde_n);
455 snew(md->tmp_v, md->mde_n);
457 snew(grpnms, md->mde_n);
459 for (i = 0; (i < md->nTC); i++)
461 ni = groups->grps[egcTC].nm_ind[i];
462 sprintf(buf, "T-%s", *(groups->grpname[ni]));
463 grpnms[i] = gmx_strdup(buf);
465 md->itemp = get_ebin_space(md->ebin, md->nTC, (const char **)grpnms,
468 if (md->etc == etcNOSEHOOVER)
470 if (md->bPrintNHChains)
472 if (md->bNHC_trotter)
474 for (i = 0; (i < md->nTC); i++)
476 ni = groups->grps[egcTC].nm_ind[i];
477 bufi = *(groups->grpname[ni]);
478 for (j = 0; (j < md->nNHC); j++)
480 sprintf(buf, "Xi-%d-%s", j, bufi);
481 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
482 sprintf(buf, "vXi-%d-%s", j, bufi);
483 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
486 md->itc = get_ebin_space(md->ebin, md->mde_n,
487 (const char **)grpnms, unit_invtime);
490 for (i = 0; (i < md->nTCP); i++)
492 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
493 for (j = 0; (j < md->nNHC); j++)
495 sprintf(buf, "Xi-%d-%s", j, bufi);
496 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
497 sprintf(buf, "vXi-%d-%s", j, bufi);
498 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
501 md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
502 (const char **)grpnms, unit_invtime);
507 for (i = 0; (i < md->nTC); i++)
509 ni = groups->grps[egcTC].nm_ind[i];
510 bufi = *(groups->grpname[ni]);
511 sprintf(buf, "Xi-%s", bufi);
512 grpnms[2*i] = gmx_strdup(buf);
513 sprintf(buf, "vXi-%s", bufi);
514 grpnms[2*i+1] = gmx_strdup(buf);
516 md->itc = get_ebin_space(md->ebin, md->mde_n,
517 (const char **)grpnms, unit_invtime);
521 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
522 md->etc == etcVRESCALE)
524 for (i = 0; (i < md->nTC); i++)
526 ni = groups->grps[egcTC].nm_ind[i];
527 sprintf(buf, "Lamb-%s", *(groups->grpname[ni]));
528 grpnms[i] = gmx_strdup(buf);
530 md->itc = get_ebin_space(md->ebin, md->mde_n, (const char **)grpnms, "");
533 for (i = 0; i < md->mde_n; i++)
539 md->nU = groups->grps[egcACC].nr;
542 snew(grpnms, 3*md->nU);
543 for (i = 0; (i < md->nU); i++)
545 ni = groups->grps[egcACC].nm_ind[i];
546 sprintf(buf, "Ux-%s", *(groups->grpname[ni]));
547 grpnms[3*i+XX] = gmx_strdup(buf);
548 sprintf(buf, "Uy-%s", *(groups->grpname[ni]));
549 grpnms[3*i+YY] = gmx_strdup(buf);
550 sprintf(buf, "Uz-%s", *(groups->grpname[ni]));
551 grpnms[3*i+ZZ] = gmx_strdup(buf);
553 md->iu = get_ebin_space(md->ebin, 3*md->nU, (const char **)grpnms, unit_vel);
559 do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
562 md->print_grpnms = nullptr;
564 /* check whether we're going to write dh histograms */
566 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
568 /* Currently dh histograms are only written with dynamics */
569 if (EI_DYNAMICS(ir->eI))
573 mde_delta_h_coll_init(md->dhc, ir);
575 md->fp_dhdl = nullptr;
576 snew(md->dE, ir->fepvals->n_lambda);
580 md->fp_dhdl = fp_dhdl;
581 snew(md->dE, ir->fepvals->n_lambda);
586 snew(md->temperatures, ir->fepvals->n_lambda);
587 for (i = 0; i < ir->fepvals->n_lambda; i++)
589 md->temperatures[i] = ir->simtempvals->temperatures[i];
595 void done_mdebin(t_mdebin *mdebin)
598 sfree(mdebin->tmp_r);
599 sfree(mdebin->tmp_v);
600 done_ebin(mdebin->ebin);
601 done_mde_delta_h_coll(mdebin->dhc);
603 sfree(mdebin->temperatures);
607 /* print a lambda vector to a string
608 fep = the inputrec's FEP input data
609 i = the index of the lambda vector
610 get_native_lambda = whether to print the native lambda
611 get_names = whether to print the names rather than the values
612 str = the pre-allocated string buffer to print to. */
613 static void print_lambda_vector(t_lambda *fep, int i,
614 gmx_bool get_native_lambda, gmx_bool get_names,
620 for (j = 0; j < efptNR; j++)
622 if (fep->separate_dvdl[j])
627 str[0] = 0; /* reset the string */
630 str += sprintf(str, "("); /* set the opening parenthesis*/
632 for (j = 0; j < efptNR; j++)
634 if (fep->separate_dvdl[j])
638 if (get_native_lambda && fep->init_lambda >= 0)
640 str += sprintf(str, "%.4f", fep->init_lambda);
644 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
649 str += sprintf(str, "%s", efpt_singular_names[j]);
651 /* print comma for the next item */
654 str += sprintf(str, ", ");
661 /* and add the closing parenthesis */
667 extern FILE *open_dhdl(const char *filename, const t_inputrec *ir,
668 const gmx_output_env_t *oenv)
671 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
672 *lambdastate = "\\lambda state";
673 char title[STRLEN], label_x[STRLEN], label_y[STRLEN];
674 int i, nps, nsets, nsets_de, nsetsbegin;
675 int n_lambda_terms = 0;
676 t_lambda *fep = ir->fepvals; /* for simplicity */
677 t_expanded *expand = ir->expandedvals;
679 char buf[STRLEN], lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
685 gmx_bool write_pV = FALSE;
687 /* count the number of different lambda terms */
688 for (i = 0; i < efptNR; i++)
690 if (fep->separate_dvdl[i])
696 if (fep->n_lambda == 0)
698 sprintf(title, "%s", dhdl);
699 sprintf(label_x, "Time (ps)");
700 sprintf(label_y, "%s (%s %s)",
701 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
705 sprintf(title, "%s and %s", dhdl, deltag);
706 sprintf(label_x, "Time (ps)");
707 sprintf(label_y, "%s and %s (%s %s)",
708 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
710 fp = gmx_fio_fopen(filename, "w+");
711 xvgr_header(fp, title, label_x, label_y, exvggtXNY, oenv);
715 bufplace = sprintf(buf, "T = %g (K) ",
718 if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
720 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
722 /* compatibility output */
723 sprintf(&(buf[bufplace]), "%s = %.4f", lambda, fep->init_lambda);
727 print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
729 print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
731 sprintf(&(buf[bufplace]), "%s %d: %s = %s",
732 lambdastate, fep->init_fep_state,
733 lambda_name_str, lambda_vec_str);
736 xvgr_subtitle(fp, buf, oenv);
740 if (fep->dhdl_derivatives == edhdlderivativesYES)
742 nsets_dhdl = n_lambda_terms;
744 /* count the number of delta_g states */
745 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
747 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
749 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
751 nsets += 1; /*add fep state for expanded ensemble */
754 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
756 nsets += 1; /* add energy to the dhdl as well */
760 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
762 nsetsextend += 1; /* for PV term, other terms possible if required for
763 the reduced potential (only needed with foreign
764 lambda, and only output when init_lambda is not
765 set in order to maintain compatibility of the
769 snew(setname, nsetsextend);
771 if (expand->elmcmove > elmcmoveNO)
773 /* state for the fep_vals, if we have alchemical sampling */
774 sprintf(buf, "%s", "Thermodynamic state");
775 setname[s] = gmx_strdup(buf);
779 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
781 switch (fep->edHdLPrintEnergy)
783 case edHdLPrintEnergyPOTENTIAL:
784 sprintf(buf, "%s (%s)", "Potential Energy", unit_energy);
786 case edHdLPrintEnergyTOTAL:
787 case edHdLPrintEnergyYES:
789 sprintf(buf, "%s (%s)", "Total Energy", unit_energy);
791 setname[s] = gmx_strdup(buf);
795 if (fep->dhdl_derivatives == edhdlderivativesYES)
797 for (i = 0; i < efptNR; i++)
799 if (fep->separate_dvdl[i])
802 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
804 /* compatibility output */
805 sprintf(buf, "%s %s %.4f", dhdl, lambda, fep->init_lambda);
809 double lam = fep->init_lambda;
810 if (fep->init_lambda < 0)
812 lam = fep->all_lambda[i][fep->init_fep_state];
814 sprintf(buf, "%s %s = %.4f", dhdl, efpt_singular_names[i],
817 setname[s] = gmx_strdup(buf);
823 if (fep->n_lambda > 0)
825 /* g_bar has to determine the lambda values used in this simulation
826 * from this xvg legend.
829 if (expand->elmcmove > elmcmoveNO)
831 nsetsbegin = 1; /* for including the expanded ensemble */
838 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
842 nsetsbegin += nsets_dhdl;
844 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
846 print_lambda_vector(fep, i, FALSE, FALSE, lambda_vec_str);
847 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
849 /* for compatible dhdl.xvg files */
850 nps = sprintf(buf, "%s %s %s", deltag, lambda, lambda_vec_str);
854 nps = sprintf(buf, "%s %s to %s", deltag, lambda, lambda_vec_str);
859 /* print the temperature for this state if doing simulated annealing */
860 sprintf(&buf[nps], "T = %g (%s)",
861 ir->simtempvals->temperatures[s-(nsetsbegin)],
864 setname[s] = gmx_strdup(buf);
869 sprintf(buf, "pV (%s)", unit_energy);
870 setname[nsetsextend-1] = gmx_strdup(buf); /* the first entry after
874 xvgr_legend(fp, nsetsextend, (const char **)setname, oenv);
876 for (s = 0; s < nsetsextend; s++)
886 static void copy_energy(t_mdebin *md, real e[], real ecpy[])
890 for (i = j = 0; (i < F_NRE); i++)
899 gmx_incons("Number of energy terms wrong");
903 void upd_mdebin(t_mdebin *md,
908 gmx_enerdata_t *enerd,
917 gmx_ekindata_t *ekind,
919 const gmx::Constraints *constr)
921 int i, j, k, kk, n, gid;
922 real crmsd[2], tmp6[6];
923 real bs[NTRICLBOXS], vol, dens, pv, enthalpy;
926 double store_dhdl[efptNR];
927 real store_energy = 0;
930 /* Do NOT use the box in the state variable, but the separate box provided
931 * as an argument. This is because we sometimes need to write the box from
932 * the last timestep to match the trajectory frames.
934 copy_energy(md, enerd->term, ecopy);
935 add_ebin(md->ebin, md->ie, md->f_nre, ecopy, bSum);
938 crmsd[0] = constr->rmsd();
939 add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
961 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
962 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
963 add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
964 add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
965 add_ebin(md->ebin, md->idens, 1, &dens, bSum);
969 /* This is pV (in kJ/mol). The pressure is the reference pressure,
970 not the instantaneous pressure */
971 pv = vol*md->ref_p/PRESFAC;
973 add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
974 enthalpy = pv + enerd->term[F_ETOT];
975 add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
980 add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
981 add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
983 add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
984 add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
985 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
986 add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
987 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
989 tmp6[0] = state->boxv[XX][XX];
990 tmp6[1] = state->boxv[YY][YY];
991 tmp6[2] = state->boxv[ZZ][ZZ];
992 tmp6[3] = state->boxv[YY][XX];
993 tmp6[4] = state->boxv[ZZ][XX];
994 tmp6[5] = state->boxv[ZZ][YY];
995 add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
999 add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1001 if (ekind && ekind->cosacc.cos_accel != 0)
1003 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1004 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1005 add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1006 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1007 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1008 *dens*gmx::square(box[ZZ][ZZ]*NANO/(2*M_PI)));
1009 add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1014 for (i = 0; (i < md->nEg); i++)
1016 for (j = i; (j < md->nEg); j++)
1018 gid = GID(i, j, md->nEg);
1019 for (k = kk = 0; (k < egNR); k++)
1023 eee[kk++] = enerd->grpp.ener[k][gid];
1026 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1034 for (i = 0; (i < md->nTC); i++)
1036 md->tmp_r[i] = ekind->tcstat[i].T;
1038 add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1040 if (md->etc == etcNOSEHOOVER)
1042 /* whether to print Nose-Hoover chains: */
1043 if (md->bPrintNHChains)
1045 if (md->bNHC_trotter)
1047 for (i = 0; (i < md->nTC); i++)
1049 for (j = 0; j < md->nNHC; j++)
1052 md->tmp_r[2*k] = state->nosehoover_xi[k];
1053 md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1056 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1060 for (i = 0; (i < md->nTCP); i++)
1062 for (j = 0; j < md->nNHC; j++)
1065 md->tmp_r[2*k] = state->nhpres_xi[k];
1066 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1069 add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1074 for (i = 0; (i < md->nTC); i++)
1076 md->tmp_r[2*i] = state->nosehoover_xi[i];
1077 md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1079 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1083 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1084 md->etc == etcVRESCALE)
1086 for (i = 0; (i < md->nTC); i++)
1088 md->tmp_r[i] = ekind->tcstat[i].lambda;
1090 add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1094 if (ekind && md->nU > 1)
1096 for (i = 0; (i < md->nU); i++)
1098 copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1100 add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1103 ebin_increase_count(md->ebin, bSum);
1105 /* BAR + thermodynamic integration values */
1106 if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1108 for (i = 0; i < enerd->n_lambda-1; i++)
1110 /* zero for simulated tempering */
1111 md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1112 if (md->temperatures != nullptr)
1114 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1115 /* is this even useful to have at all? */
1116 md->dE[i] += (md->temperatures[i]/
1117 md->temperatures[state->fep_state]-1.0)*
1118 enerd->term[F_EKIN];
1124 fprintf(md->fp_dhdl, "%.4f", time);
1125 /* the current free energy state */
1127 /* print the current state if we are doing expanded ensemble */
1128 if (expand->elmcmove > elmcmoveNO)
1130 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1132 /* total energy (for if the temperature changes */
1134 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1136 switch (fep->edHdLPrintEnergy)
1138 case edHdLPrintEnergyPOTENTIAL:
1139 store_energy = enerd->term[F_EPOT];
1141 case edHdLPrintEnergyTOTAL:
1142 case edHdLPrintEnergyYES:
1144 store_energy = enerd->term[F_ETOT];
1146 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1149 if (fep->dhdl_derivatives == edhdlderivativesYES)
1151 for (i = 0; i < efptNR; i++)
1153 if (fep->separate_dvdl[i])
1155 /* assumes F_DVDL is first */
1156 fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1160 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1162 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1166 (md->epc != epcNO) &&
1167 (enerd->n_lambda > 0) &&
1168 (fep->init_lambda < 0))
1170 fprintf(md->fp_dhdl, " %#.8g", pv); /* PV term only needed when
1171 there are alternate state
1172 lambda and we're not in
1173 compatibility mode */
1175 fprintf(md->fp_dhdl, "\n");
1176 /* and the binary free energy output */
1178 if (md->dhc && bDoDHDL)
1181 for (i = 0; i < efptNR; i++)
1183 if (fep->separate_dvdl[i])
1185 /* assumes F_DVDL is first */
1186 store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1190 store_energy = enerd->term[F_ETOT];
1191 /* store_dh is dE */
1192 mde_delta_h_coll_add_dh(md->dhc,
1193 (double)state->fep_state,
1197 md->dE + fep->lambda_start_n,
1204 void upd_mdebin_step(t_mdebin *md)
1206 ebin_increase_count(md->ebin, FALSE);
1209 static void npr(FILE *log, int n, char c)
1211 for (; (n > 0); n--)
1213 fprintf(log, "%c", c);
1217 static void pprint(FILE *log, const char *s, t_mdebin *md)
1221 char buf1[22], buf2[22];
1224 fprintf(log, "\t<====== ");
1225 npr(log, slen, CHAR);
1226 fprintf(log, " ==>\n");
1227 fprintf(log, "\t<==== %s ====>\n", s);
1228 fprintf(log, "\t<== ");
1229 npr(log, slen, CHAR);
1230 fprintf(log, " ======>\n\n");
1232 fprintf(log, "\tStatistics over %s steps using %s frames\n",
1233 gmx_step_str(md->ebin->nsteps_sim, buf1),
1234 gmx_step_str(md->ebin->nsum_sim, buf2));
1238 void print_ebin_header(FILE *log, gmx_int64_t steps, double time)
1242 fprintf(log, " %12s %12s\n"
1244 "Step", "Time", gmx_step_str(steps, buf), time);
1247 // TODO It is too many responsibilities for this function to handle
1248 // both .edr and .log output for both per-time and time-average data.
1249 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1251 gmx_int64_t step, double time,
1253 t_mdebin *md, t_fcdata *fcd,
1254 gmx_groups_t *groups, t_grpopts *opts,
1257 /*static char **grpnms=NULL;*/
1259 int i, j, n, ni, nj, b;
1261 real *disre_rm3tav, *disre_rt;
1263 /* these are for the old-style blocks (1 subblock, only reals), because
1264 there can be only one per ID for these */
1271 if (mode == eprAVER && md->ebin->nsum_sim <= 0)
1275 fprintf(log, "Not enough data recorded to report energy averages\n");
1286 fr.nsteps = md->ebin->nsteps;
1287 fr.dt = md->delta_t;
1288 fr.nsum = md->ebin->nsum;
1289 fr.nre = (bEne) ? md->ebin->nener : 0;
1290 fr.ener = md->ebin->e;
1291 ndisre = bDR ? fcd->disres.npair : 0;
1292 disre_rm3tav = fcd->disres.rm3tav;
1293 disre_rt = fcd->disres.rt;
1294 /* Optional additional old-style (real-only) blocks. */
1295 for (i = 0; i < enxNR; i++)
1299 if (fcd->orires.nr > 0 && bOR)
1301 diagonalize_orires_tensors(&(fcd->orires));
1302 nr[enxOR] = fcd->orires.nr;
1303 block[enxOR] = fcd->orires.otav;
1305 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1307 block[enxORI] = fcd->orires.oinsl;
1308 id[enxORI] = enxORI;
1309 nr[enxORT] = fcd->orires.nex*12;
1310 block[enxORT] = fcd->orires.eig;
1311 id[enxORT] = enxORT;
1314 /* whether we are going to wrte anything out: */
1315 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1318 /* the old-style blocks go first */
1320 for (i = 0; i < enxNR; i++)
1327 add_blocks_enxframe(&fr, fr.nblock);
1328 for (b = 0; b < fr.nblock; b++)
1330 add_subblocks_enxblock(&(fr.block[b]), 1);
1331 fr.block[b].id = id[b];
1332 fr.block[b].sub[0].nr = nr[b];
1334 fr.block[b].sub[0].type = xdr_datatype_float;
1335 fr.block[b].sub[0].fval = block[b];
1337 fr.block[b].sub[0].type = xdr_datatype_double;
1338 fr.block[b].sub[0].dval = block[b];
1342 /* check for disre block & fill it. */
1347 add_blocks_enxframe(&fr, fr.nblock);
1349 add_subblocks_enxblock(&(fr.block[db]), 2);
1350 fr.block[db].id = enxDISRE;
1351 fr.block[db].sub[0].nr = ndisre;
1352 fr.block[db].sub[1].nr = ndisre;
1354 fr.block[db].sub[0].type = xdr_datatype_float;
1355 fr.block[db].sub[1].type = xdr_datatype_float;
1356 fr.block[db].sub[0].fval = disre_rt;
1357 fr.block[db].sub[1].fval = disre_rm3tav;
1359 fr.block[db].sub[0].type = xdr_datatype_double;
1360 fr.block[db].sub[1].type = xdr_datatype_double;
1361 fr.block[db].sub[0].dval = disre_rt;
1362 fr.block[db].sub[1].dval = disre_rm3tav;
1365 /* here we can put new-style blocks */
1367 /* Free energy perturbation blocks */
1370 mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1373 /* we can now free & reset the data in the blocks */
1376 mde_delta_h_coll_reset(md->dhc);
1379 /* AWH bias blocks. */
1380 if (awh != nullptr) // TODO: add boolean in t_mdebin. See in mdebin.h.
1382 awh->writeToEnergyFrame(step, &fr);
1385 /* do the actual I/O */
1386 do_enx(fp_ene, &fr);
1389 /* We have stored the sums, so reset the sum history */
1390 reset_ebin_sums(md->ebin);
1398 pprint(log, "A V E R A G E S", md);
1404 pprint(log, "R M S - F L U C T U A T I O N S", md);
1408 gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1413 for (i = 0; i < opts->ngtc; i++)
1415 if (opts->annealing[i] != eannNO)
1417 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1418 *(groups->grpname[groups->grps[egcTC].nm_ind[i]]),
1422 if (mode == eprNORMAL && fcd->orires.nr > 0)
1424 print_orires_log(log, &(fcd->orires));
1426 fprintf(log, " Energies (%s)\n", unit_energy);
1427 pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1430 if (mode == eprAVER)
1434 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1440 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1441 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1443 fprintf(log, " Force Virial (%s)\n", unit_energy);
1444 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1447 fprintf(log, " Total Virial (%s)\n", unit_energy);
1448 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1450 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1451 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1455 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1456 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1462 if (md->print_grpnms == nullptr)
1464 snew(md->print_grpnms, md->nE);
1466 for (i = 0; (i < md->nEg); i++)
1468 ni = groups->grps[egcENER].nm_ind[i];
1469 for (j = i; (j < md->nEg); j++)
1471 nj = groups->grps[egcENER].nm_ind[j];
1472 sprintf(buf, "%s-%s", *(groups->grpname[ni]),
1473 *(groups->grpname[nj]));
1474 md->print_grpnms[n++] = gmx_strdup(buf);
1478 sprintf(buf, "Epot (%s)", unit_energy);
1479 fprintf(log, "%15s ", buf);
1480 for (i = 0; (i < egNR); i++)
1484 fprintf(log, "%12s ", egrp_nm[i]);
1488 for (i = 0; (i < md->nE); i++)
1490 fprintf(log, "%15s", md->print_grpnms[i]);
1491 pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1498 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1503 fprintf(log, "%15s %12s %12s %12s\n",
1504 "Group", "Ux", "Uy", "Uz");
1505 for (i = 0; (i < md->nU); i++)
1507 ni = groups->grps[egcACC].nm_ind[i];
1508 fprintf(log, "%15s", *groups->grpname[ni]);
1509 pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1518 void update_energyhistory(energyhistory_t * enerhist, const t_mdebin * mdebin)
1520 const t_ebin * const ebin = mdebin->ebin;
1522 enerhist->nsteps = ebin->nsteps;
1523 enerhist->nsum = ebin->nsum;
1524 enerhist->nsteps_sim = ebin->nsteps_sim;
1525 enerhist->nsum_sim = ebin->nsum_sim;
1529 /* This will only actually resize the first time */
1530 enerhist->ener_ave.resize(ebin->nener);
1531 enerhist->ener_sum.resize(ebin->nener);
1533 for (int i = 0; i < ebin->nener; i++)
1535 enerhist->ener_ave[i] = ebin->e[i].eav;
1536 enerhist->ener_sum[i] = ebin->e[i].esum;
1540 if (ebin->nsum_sim > 0)
1542 /* This will only actually resize the first time */
1543 enerhist->ener_sum_sim.resize(ebin->nener);
1545 for (int i = 0; i < ebin->nener; i++)
1547 enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1552 mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1556 void restore_energyhistory_from_state(t_mdebin * mdebin,
1557 const energyhistory_t * enerhist)
1559 unsigned int nener = static_cast<unsigned int>(mdebin->ebin->nener);
1561 if ((enerhist->nsum > 0 && nener != enerhist->ener_sum.size()) ||
1562 (enerhist->nsum_sim > 0 && nener != enerhist->ener_sum_sim.size()))
1564 gmx_fatal(FARGS, "Mismatch between number of energies in run input (%d) and checkpoint file (%u or %u).",
1565 nener, enerhist->ener_sum.size(), enerhist->ener_sum_sim.size());
1568 mdebin->ebin->nsteps = enerhist->nsteps;
1569 mdebin->ebin->nsum = enerhist->nsum;
1570 mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1571 mdebin->ebin->nsum_sim = enerhist->nsum_sim;
1573 for (int i = 0; i < mdebin->ebin->nener; i++)
1575 mdebin->ebin->e[i].eav =
1576 (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1577 mdebin->ebin->e[i].esum =
1578 (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1579 mdebin->ebin->e_sim[i].esum =
1580 (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1584 mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist->deltaHForeignLambdas.get());