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 const char *egrp_nm[egNR+1] = {
98 "Coul-SR", "LJ-SR", "Buck-SR",
99 "Coul-14", "LJ-14", nullptr
102 t_mdebin *init_mdebin(ener_file_t fp_ene,
103 const gmx_mtop_t *mtop,
104 const t_inputrec *ir,
107 const char *ener_nm[F_NRE];
108 static const char *vir_nm[] = {
109 "Vir-XX", "Vir-XY", "Vir-XZ",
110 "Vir-YX", "Vir-YY", "Vir-YZ",
111 "Vir-ZX", "Vir-ZY", "Vir-ZZ"
113 static const char *sv_nm[] = {
114 "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
115 "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
116 "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ"
118 static const char *fv_nm[] = {
119 "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
120 "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
121 "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ"
123 static const char *pres_nm[] = {
124 "Pres-XX", "Pres-XY", "Pres-XZ",
125 "Pres-YX", "Pres-YY", "Pres-YZ",
126 "Pres-ZX", "Pres-ZY", "Pres-ZZ"
128 static const char *surft_nm[] = {
131 static const char *mu_nm[] = {
132 "Mu-X", "Mu-Y", "Mu-Z"
134 static const char *vcos_nm[] = {
137 static const char *visc_nm[] = {
140 static const char *baro_nm[] = {
144 const gmx_groups_t *groups;
149 int i, j, ni, nj, n, k, kk, ncon, nset;
154 if (EI_DYNAMICS(ir->eI))
156 md->delta_t = ir->delta_t;
163 groups = &mtop->groups;
165 bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
166 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
167 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
169 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
170 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
171 md->bConstr = (ncon > 0 || nset > 0);
172 md->bConstrVir = FALSE;
175 if (ncon > 0 && ir->eConstrAlg == econtLINCS)
179 md->bConstrVir = (getenv("GMX_CONSTRAINTVIR") != nullptr);
186 /* Energy monitoring */
187 for (i = 0; i < egNR; i++)
189 md->bEInd[i] = FALSE;
192 for (i = 0; i < F_NRE; i++)
194 md->bEner[i] = FALSE;
197 md->bEner[i] = !bBHAM;
199 else if (i == F_BHAM)
201 md->bEner[i] = bBHAM;
205 md->bEner[i] = ir->bQMMM;
207 else if (i == F_RF_EXCL)
209 md->bEner[i] = (EEL_RF(ir->coulombtype) && ir->cutoff_scheme == ecutsGROUP);
211 else if (i == F_COUL_RECIP)
213 md->bEner[i] = EEL_FULL(ir->coulombtype);
215 else if (i == F_LJ_RECIP)
217 md->bEner[i] = EVDW_PME(ir->vdwtype);
219 else if (i == F_LJ14)
223 else if (i == F_COUL14)
227 else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
229 md->bEner[i] = FALSE;
231 else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
232 (i == F_DVDL_VDW && ir->fepvals->separate_dvdl[efptVDW]) ||
233 (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
234 (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
235 (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
236 (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
238 md->bEner[i] = (ir->efep != efepNO);
240 else if ((interaction_function[i].flags & IF_VSITE) ||
241 (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
243 md->bEner[i] = FALSE;
245 else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES) || (i == F_EQM))
249 else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
251 md->bEner[i] = EI_DYNAMICS(ir->eI);
253 else if (i == F_DISPCORR || i == F_PDISPCORR)
255 md->bEner[i] = (ir->eDispCorr != edispcNO);
257 else if (i == F_DISRESVIOL)
259 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
261 else if (i == F_ORIRESDEV)
263 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
265 else if (i == F_CONNBONDS)
267 md->bEner[i] = FALSE;
269 else if (i == F_COM_PULL)
271 md->bEner[i] = ((ir->bPull && pull_have_potential(ir->pull_work)) ||
274 else if (i == F_ECONSERVED)
276 md->bEner[i] = (integratorHasConservedEnergyQuantity(ir));
280 md->bEner[i] = (gmx_mtop_ftype_count(mtop, i) > 0);
285 for (i = 0; i < F_NRE; i++)
289 ener_nm[md->f_nre] = interaction_function[i].longname;
295 md->bDiagPres = !TRICLINIC(ir->ref_p);
296 md->ref_p = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
297 md->bTricl = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
298 md->bDynBox = inputrecDynamicBox(ir);
300 md->bNHC_trotter = inputrecNvtTrotter(ir);
301 md->bPrintNHChains = ir->bPrintNHChains;
302 md->bMTTK = (inputrecNptTrotter(ir) || inputrecNphTrotter(ir));
303 md->bMu = inputrecNeedMutot(ir);
305 md->ebin = mk_ebin();
306 /* Pass NULL for unit to let get_ebin_space determine the units
307 * for interaction_function[i].longname
309 md->ie = get_ebin_space(md->ebin, md->f_nre, ener_nm, nullptr);
312 /* This should be called directly after the call for md->ie,
313 * such that md->iconrmsd follows directly in the list.
315 md->iconrmsd = get_ebin_space(md->ebin, md->nCrmsd, conrmsd_nm, "");
319 md->ib = get_ebin_space(md->ebin,
320 md->bTricl ? NTRICLBOXS : NBOXS,
321 md->bTricl ? tricl_boxs_nm : boxs_nm,
323 md->ivol = get_ebin_space(md->ebin, 1, vol_nm, unit_volume);
324 md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
327 md->ipv = get_ebin_space(md->ebin, 1, pv_nm, unit_energy);
328 md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm, unit_energy);
333 md->isvir = get_ebin_space(md->ebin, asize(sv_nm), sv_nm, unit_energy);
334 md->ifvir = get_ebin_space(md->ebin, asize(fv_nm), fv_nm, unit_energy);
336 md->ivir = get_ebin_space(md->ebin, asize(vir_nm), vir_nm, unit_energy);
337 md->ipres = get_ebin_space(md->ebin, asize(pres_nm), pres_nm, unit_pres_bar);
338 md->isurft = get_ebin_space(md->ebin, asize(surft_nm), surft_nm,
340 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
342 md->ipc = get_ebin_space(md->ebin, md->bTricl ? 6 : 3,
343 boxvel_nm, unit_vel);
347 md->imu = get_ebin_space(md->ebin, asize(mu_nm), mu_nm, unit_dipole_D);
349 if (ir->cos_accel != 0)
351 md->ivcos = get_ebin_space(md->ebin, asize(vcos_nm), vcos_nm, unit_vel);
352 md->ivisc = get_ebin_space(md->ebin, asize(visc_nm), visc_nm,
356 /* Energy monitoring */
357 for (i = 0; i < egNR; i++)
359 md->bEInd[i] = FALSE;
361 md->bEInd[egCOULSR] = TRUE;
362 md->bEInd[egLJSR ] = TRUE;
366 md->bEInd[egLJSR] = FALSE;
367 md->bEInd[egBHAMSR] = TRUE;
371 md->bEInd[egLJ14] = TRUE;
372 md->bEInd[egCOUL14] = TRUE;
375 for (i = 0; (i < egNR); i++)
383 n = groups->grps[egcENER].nr;
385 md->nE = (n*(n+1))/2;
387 snew(md->igrp, md->nE);
392 for (k = 0; (k < md->nEc); k++)
394 snew(gnm[k], STRLEN);
396 for (i = 0; (i < groups->grps[egcENER].nr); i++)
398 ni = groups->grps[egcENER].nm_ind[i];
399 for (j = i; (j < groups->grps[egcENER].nr); j++)
401 nj = groups->grps[egcENER].nm_ind[j];
402 for (k = kk = 0; (k < egNR); k++)
406 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
407 *(groups->grpname[ni]), *(groups->grpname[nj]));
411 md->igrp[n] = get_ebin_space(md->ebin, md->nEc,
416 for (k = 0; (k < md->nEc); k++)
424 gmx_incons("Number of energy terms wrong");
428 md->nTC = groups->grps[egcTC].nr;
429 md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
432 md->nTCP = 1; /* assume only one possible coupling system for barostat
439 if (md->etc == etcNOSEHOOVER)
441 if (md->bNHC_trotter)
443 md->mde_n = 2*md->nNHC*md->nTC;
447 md->mde_n = 2*md->nTC;
449 if (md->epc == epcMTTK)
451 md->mdeb_n = 2*md->nNHC*md->nTCP;
460 snew(md->tmp_r, md->mde_n);
461 snew(md->tmp_v, md->mde_n);
463 snew(grpnms, md->mde_n);
465 for (i = 0; (i < md->nTC); i++)
467 ni = groups->grps[egcTC].nm_ind[i];
468 sprintf(buf, "T-%s", *(groups->grpname[ni]));
469 grpnms[i] = gmx_strdup(buf);
471 md->itemp = get_ebin_space(md->ebin, md->nTC, grpnms,
474 if (md->etc == etcNOSEHOOVER)
476 if (md->bPrintNHChains)
478 if (md->bNHC_trotter)
480 for (i = 0; (i < md->nTC); i++)
482 ni = groups->grps[egcTC].nm_ind[i];
483 bufi = *(groups->grpname[ni]);
484 for (j = 0; (j < md->nNHC); j++)
486 sprintf(buf, "Xi-%d-%s", j, bufi);
487 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
488 sprintf(buf, "vXi-%d-%s", j, bufi);
489 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
492 md->itc = get_ebin_space(md->ebin, md->mde_n,
493 grpnms, unit_invtime);
496 for (i = 0; (i < md->nTCP); i++)
498 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
499 for (j = 0; (j < md->nNHC); j++)
501 sprintf(buf, "Xi-%d-%s", j, bufi);
502 grpnms[2*(i*md->nNHC+j)] = gmx_strdup(buf);
503 sprintf(buf, "vXi-%d-%s", j, bufi);
504 grpnms[2*(i*md->nNHC+j)+1] = gmx_strdup(buf);
507 md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
508 grpnms, unit_invtime);
513 for (i = 0; (i < md->nTC); i++)
515 ni = groups->grps[egcTC].nm_ind[i];
516 bufi = *(groups->grpname[ni]);
517 sprintf(buf, "Xi-%s", bufi);
518 grpnms[2*i] = gmx_strdup(buf);
519 sprintf(buf, "vXi-%s", bufi);
520 grpnms[2*i+1] = gmx_strdup(buf);
522 md->itc = get_ebin_space(md->ebin, md->mde_n,
523 grpnms, unit_invtime);
527 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
528 md->etc == etcVRESCALE)
530 for (i = 0; (i < md->nTC); i++)
532 ni = groups->grps[egcTC].nm_ind[i];
533 sprintf(buf, "Lamb-%s", *(groups->grpname[ni]));
534 grpnms[i] = gmx_strdup(buf);
536 md->itc = get_ebin_space(md->ebin, md->mde_n, grpnms, "");
539 for (i = 0; i < md->mde_n; i++)
545 md->nU = groups->grps[egcACC].nr;
548 snew(grpnms, 3*md->nU);
549 for (i = 0; (i < md->nU); i++)
551 ni = groups->grps[egcACC].nm_ind[i];
552 sprintf(buf, "Ux-%s", *(groups->grpname[ni]));
553 grpnms[3*i+XX] = gmx_strdup(buf);
554 sprintf(buf, "Uy-%s", *(groups->grpname[ni]));
555 grpnms[3*i+YY] = gmx_strdup(buf);
556 sprintf(buf, "Uz-%s", *(groups->grpname[ni]));
557 grpnms[3*i+ZZ] = gmx_strdup(buf);
559 md->iu = get_ebin_space(md->ebin, 3*md->nU, grpnms, unit_vel);
565 do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
568 md->print_grpnms = nullptr;
570 /* check whether we're going to write dh histograms */
572 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
574 /* Currently dh histograms are only written with dynamics */
575 if (EI_DYNAMICS(ir->eI))
579 mde_delta_h_coll_init(md->dhc, ir);
581 md->fp_dhdl = nullptr;
582 snew(md->dE, ir->fepvals->n_lambda);
586 md->fp_dhdl = fp_dhdl;
587 snew(md->dE, ir->fepvals->n_lambda);
592 snew(md->temperatures, ir->fepvals->n_lambda);
593 for (i = 0; i < ir->fepvals->n_lambda; i++)
595 md->temperatures[i] = ir->simtempvals->temperatures[i];
601 void done_mdebin(t_mdebin *mdebin)
604 sfree(mdebin->tmp_r);
605 sfree(mdebin->tmp_v);
606 done_ebin(mdebin->ebin);
607 done_mde_delta_h_coll(mdebin->dhc);
609 sfree(mdebin->temperatures);
613 /* print a lambda vector to a string
614 fep = the inputrec's FEP input data
615 i = the index of the lambda vector
616 get_native_lambda = whether to print the native lambda
617 get_names = whether to print the names rather than the values
618 str = the pre-allocated string buffer to print to. */
619 static void print_lambda_vector(t_lambda *fep, int i,
620 gmx_bool get_native_lambda, gmx_bool get_names,
626 for (j = 0; j < efptNR; j++)
628 if (fep->separate_dvdl[j])
633 str[0] = 0; /* reset the string */
636 str += sprintf(str, "("); /* set the opening parenthesis*/
638 for (j = 0; j < efptNR; j++)
640 if (fep->separate_dvdl[j])
644 if (get_native_lambda && fep->init_lambda >= 0)
646 str += sprintf(str, "%.4f", fep->init_lambda);
650 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
655 str += sprintf(str, "%s", efpt_singular_names[j]);
657 /* print comma for the next item */
660 str += sprintf(str, ", ");
667 /* and add the closing parenthesis */
673 extern FILE *open_dhdl(const char *filename, const t_inputrec *ir,
674 const gmx_output_env_t *oenv)
677 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
678 *lambdastate = "\\lambda state";
679 char title[STRLEN], label_x[STRLEN], label_y[STRLEN];
680 int i, nps, nsets, nsets_de, nsetsbegin;
681 int n_lambda_terms = 0;
682 t_lambda *fep = ir->fepvals; /* for simplicity */
683 t_expanded *expand = ir->expandedvals;
685 char buf[STRLEN], lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
691 gmx_bool write_pV = FALSE;
693 /* count the number of different lambda terms */
694 for (i = 0; i < efptNR; i++)
696 if (fep->separate_dvdl[i])
702 if (fep->n_lambda == 0)
704 sprintf(title, "%s", dhdl);
705 sprintf(label_x, "Time (ps)");
706 sprintf(label_y, "%s (%s %s)",
707 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
711 sprintf(title, "%s and %s", dhdl, deltag);
712 sprintf(label_x, "Time (ps)");
713 sprintf(label_y, "%s and %s (%s %s)",
714 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
716 fp = gmx_fio_fopen(filename, "w+");
717 xvgr_header(fp, title, label_x, label_y, exvggtXNY, oenv);
721 bufplace = sprintf(buf, "T = %g (K) ",
724 if ((ir->efep != efepSLOWGROWTH) && (ir->efep != efepEXPANDED))
726 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
728 /* compatibility output */
729 sprintf(&(buf[bufplace]), "%s = %.4f", lambda, fep->init_lambda);
733 print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
735 print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
737 sprintf(&(buf[bufplace]), "%s %d: %s = %s",
738 lambdastate, fep->init_fep_state,
739 lambda_name_str, lambda_vec_str);
742 xvgr_subtitle(fp, buf, oenv);
746 if (fep->dhdl_derivatives == edhdlderivativesYES)
748 nsets_dhdl = n_lambda_terms;
750 /* count the number of delta_g states */
751 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
753 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
755 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
757 nsets += 1; /*add fep state for expanded ensemble */
760 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
762 nsets += 1; /* add energy to the dhdl as well */
766 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
768 nsetsextend += 1; /* for PV term, other terms possible if required for
769 the reduced potential (only needed with foreign
770 lambda, and only output when init_lambda is not
771 set in order to maintain compatibility of the
775 snew(setname, nsetsextend);
777 if (expand->elmcmove > elmcmoveNO)
779 /* state for the fep_vals, if we have alchemical sampling */
780 sprintf(buf, "%s", "Thermodynamic state");
781 setname[s] = gmx_strdup(buf);
785 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
787 switch (fep->edHdLPrintEnergy)
789 case edHdLPrintEnergyPOTENTIAL:
790 sprintf(buf, "%s (%s)", "Potential Energy", unit_energy);
792 case edHdLPrintEnergyTOTAL:
793 case edHdLPrintEnergyYES:
795 sprintf(buf, "%s (%s)", "Total Energy", unit_energy);
797 setname[s] = gmx_strdup(buf);
801 if (fep->dhdl_derivatives == edhdlderivativesYES)
803 for (i = 0; i < efptNR; i++)
805 if (fep->separate_dvdl[i])
808 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
810 /* compatibility output */
811 sprintf(buf, "%s %s %.4f", dhdl, lambda, fep->init_lambda);
815 double lam = fep->init_lambda;
816 if (fep->init_lambda < 0)
818 lam = fep->all_lambda[i][fep->init_fep_state];
820 sprintf(buf, "%s %s = %.4f", dhdl, efpt_singular_names[i],
823 setname[s] = gmx_strdup(buf);
829 if (fep->n_lambda > 0)
831 /* g_bar has to determine the lambda values used in this simulation
832 * from this xvg legend.
835 if (expand->elmcmove > elmcmoveNO)
837 nsetsbegin = 1; /* for including the expanded ensemble */
844 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
848 nsetsbegin += nsets_dhdl;
850 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
852 print_lambda_vector(fep, i, FALSE, FALSE, lambda_vec_str);
853 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
855 /* for compatible dhdl.xvg files */
856 nps = sprintf(buf, "%s %s %s", deltag, lambda, lambda_vec_str);
860 nps = sprintf(buf, "%s %s to %s", deltag, lambda, lambda_vec_str);
865 /* print the temperature for this state if doing simulated annealing */
866 sprintf(&buf[nps], "T = %g (%s)",
867 ir->simtempvals->temperatures[s-(nsetsbegin)],
870 setname[s] = gmx_strdup(buf);
875 sprintf(buf, "pV (%s)", unit_energy);
876 setname[nsetsextend-1] = gmx_strdup(buf); /* the first entry after
880 xvgr_legend(fp, nsetsextend, setname, oenv);
882 for (s = 0; s < nsetsextend; s++)
892 static void copy_energy(t_mdebin *md, const real e[], real ecpy[])
896 for (i = j = 0; (i < F_NRE); i++)
905 gmx_incons("Number of energy terms wrong");
909 void upd_mdebin(t_mdebin *md,
914 gmx_enerdata_t *enerd,
923 gmx_ekindata_t *ekind,
925 const gmx::Constraints *constr)
927 int i, j, k, kk, n, gid;
928 real crmsd[2], tmp6[6];
929 real bs[NTRICLBOXS], vol, dens, pv, enthalpy;
932 double store_dhdl[efptNR];
933 real store_energy = 0;
936 /* Do NOT use the box in the state variable, but the separate box provided
937 * as an argument. This is because we sometimes need to write the box from
938 * the last timestep to match the trajectory frames.
940 copy_energy(md, enerd->term, ecopy);
941 add_ebin(md->ebin, md->ie, md->f_nre, ecopy, bSum);
944 crmsd[0] = constr->rmsd();
945 add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
967 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
968 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
969 add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
970 add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
971 add_ebin(md->ebin, md->idens, 1, &dens, bSum);
975 /* This is pV (in kJ/mol). The pressure is the reference pressure,
976 not the instantaneous pressure */
977 pv = vol*md->ref_p/PRESFAC;
979 add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
980 enthalpy = pv + enerd->term[F_ETOT];
981 add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
986 add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
987 add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
989 add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
990 add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
991 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
992 add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
993 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
995 tmp6[0] = state->boxv[XX][XX];
996 tmp6[1] = state->boxv[YY][YY];
997 tmp6[2] = state->boxv[ZZ][ZZ];
998 tmp6[3] = state->boxv[YY][XX];
999 tmp6[4] = state->boxv[ZZ][XX];
1000 tmp6[5] = state->boxv[ZZ][YY];
1001 add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1005 add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1007 if (ekind && ekind->cosacc.cos_accel != 0)
1009 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1010 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1011 add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1012 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1013 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1014 *dens*gmx::square(box[ZZ][ZZ]*NANO/(2*M_PI)));
1015 add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1020 for (i = 0; (i < md->nEg); i++)
1022 for (j = i; (j < md->nEg); j++)
1024 gid = GID(i, j, md->nEg);
1025 for (k = kk = 0; (k < egNR); k++)
1029 eee[kk++] = enerd->grpp.ener[k][gid];
1032 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1040 for (i = 0; (i < md->nTC); i++)
1042 md->tmp_r[i] = ekind->tcstat[i].T;
1044 add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1046 if (md->etc == etcNOSEHOOVER)
1048 /* whether to print Nose-Hoover chains: */
1049 if (md->bPrintNHChains)
1051 if (md->bNHC_trotter)
1053 for (i = 0; (i < md->nTC); i++)
1055 for (j = 0; j < md->nNHC; j++)
1058 md->tmp_r[2*k] = state->nosehoover_xi[k];
1059 md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1062 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1066 for (i = 0; (i < md->nTCP); i++)
1068 for (j = 0; j < md->nNHC; j++)
1071 md->tmp_r[2*k] = state->nhpres_xi[k];
1072 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1075 add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1080 for (i = 0; (i < md->nTC); i++)
1082 md->tmp_r[2*i] = state->nosehoover_xi[i];
1083 md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1085 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1089 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1090 md->etc == etcVRESCALE)
1092 for (i = 0; (i < md->nTC); i++)
1094 md->tmp_r[i] = ekind->tcstat[i].lambda;
1096 add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1100 if (ekind && md->nU > 1)
1102 for (i = 0; (i < md->nU); i++)
1104 copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1106 add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1109 ebin_increase_count(md->ebin, bSum);
1111 /* BAR + thermodynamic integration values */
1112 if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1114 for (i = 0; i < enerd->n_lambda-1; i++)
1116 /* zero for simulated tempering */
1117 md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1118 if (md->temperatures != nullptr)
1120 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1121 /* is this even useful to have at all? */
1122 md->dE[i] += (md->temperatures[i]/
1123 md->temperatures[state->fep_state]-1.0)*
1124 enerd->term[F_EKIN];
1130 fprintf(md->fp_dhdl, "%.4f", time);
1131 /* the current free energy state */
1133 /* print the current state if we are doing expanded ensemble */
1134 if (expand->elmcmove > elmcmoveNO)
1136 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1138 /* total energy (for if the temperature changes */
1140 if (fep->edHdLPrintEnergy != edHdLPrintEnergyNO)
1142 switch (fep->edHdLPrintEnergy)
1144 case edHdLPrintEnergyPOTENTIAL:
1145 store_energy = enerd->term[F_EPOT];
1147 case edHdLPrintEnergyTOTAL:
1148 case edHdLPrintEnergyYES:
1150 store_energy = enerd->term[F_ETOT];
1152 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1155 if (fep->dhdl_derivatives == edhdlderivativesYES)
1157 for (i = 0; i < efptNR; i++)
1159 if (fep->separate_dvdl[i])
1161 /* assumes F_DVDL is first */
1162 fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1166 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1168 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1172 (md->epc != epcNO) &&
1173 (enerd->n_lambda > 0) &&
1174 (fep->init_lambda < 0))
1176 fprintf(md->fp_dhdl, " %#.8g", pv); /* PV term only needed when
1177 there are alternate state
1178 lambda and we're not in
1179 compatibility mode */
1181 fprintf(md->fp_dhdl, "\n");
1182 /* and the binary free energy output */
1184 if (md->dhc && bDoDHDL)
1187 for (i = 0; i < efptNR; i++)
1189 if (fep->separate_dvdl[i])
1191 /* assumes F_DVDL is first */
1192 store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1196 store_energy = enerd->term[F_ETOT];
1197 /* store_dh is dE */
1198 mde_delta_h_coll_add_dh(md->dhc,
1199 static_cast<double>(state->fep_state),
1203 md->dE + fep->lambda_start_n,
1210 void upd_mdebin_step(t_mdebin *md)
1212 ebin_increase_count(md->ebin, FALSE);
1215 static void npr(FILE *log, int n, char c)
1217 for (; (n > 0); n--)
1219 fprintf(log, "%c", c);
1223 static void pprint(FILE *log, const char *s, t_mdebin *md)
1227 char buf1[22], buf2[22];
1230 fprintf(log, "\t<====== ");
1231 npr(log, slen, CHAR);
1232 fprintf(log, " ==>\n");
1233 fprintf(log, "\t<==== %s ====>\n", s);
1234 fprintf(log, "\t<== ");
1235 npr(log, slen, CHAR);
1236 fprintf(log, " ======>\n\n");
1238 fprintf(log, "\tStatistics over %s steps using %s frames\n",
1239 gmx_step_str(md->ebin->nsteps_sim, buf1),
1240 gmx_step_str(md->ebin->nsum_sim, buf2));
1244 void print_ebin_header(FILE *log, int64_t steps, double time)
1248 fprintf(log, " %12s %12s\n"
1250 "Step", "Time", gmx_step_str(steps, buf), time);
1253 // TODO It is too many responsibilities for this function to handle
1254 // both .edr and .log output for both per-time and time-average data.
1255 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1257 int64_t step, double time,
1259 t_mdebin *md, t_fcdata *fcd,
1260 gmx_groups_t *groups, t_grpopts *opts,
1263 /*static char **grpnms=NULL;*/
1265 int i, j, n, ni, nj, b;
1267 real *disre_rm3tav, *disre_rt;
1269 /* these are for the old-style blocks (1 subblock, only reals), because
1270 there can be only one per ID for these */
1277 if (mode == eprAVER && md->ebin->nsum_sim <= 0)
1281 fprintf(log, "Not enough data recorded to report energy averages\n");
1292 fr.nsteps = md->ebin->nsteps;
1293 fr.dt = md->delta_t;
1294 fr.nsum = md->ebin->nsum;
1295 fr.nre = (bEne) ? md->ebin->nener : 0;
1296 fr.ener = md->ebin->e;
1297 ndisre = bDR ? fcd->disres.npair : 0;
1298 disre_rm3tav = fcd->disres.rm3tav;
1299 disre_rt = fcd->disres.rt;
1300 /* Optional additional old-style (real-only) blocks. */
1301 for (i = 0; i < enxNR; i++)
1305 if (fcd->orires.nr > 0 && bOR)
1307 diagonalize_orires_tensors(&(fcd->orires));
1308 nr[enxOR] = fcd->orires.nr;
1309 block[enxOR] = fcd->orires.otav;
1311 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1313 block[enxORI] = fcd->orires.oinsl;
1314 id[enxORI] = enxORI;
1315 nr[enxORT] = fcd->orires.nex*12;
1316 block[enxORT] = fcd->orires.eig;
1317 id[enxORT] = enxORT;
1320 /* whether we are going to wrte anything out: */
1321 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1324 /* the old-style blocks go first */
1326 for (i = 0; i < enxNR; i++)
1333 add_blocks_enxframe(&fr, fr.nblock);
1334 for (b = 0; b < fr.nblock; b++)
1336 add_subblocks_enxblock(&(fr.block[b]), 1);
1337 fr.block[b].id = id[b];
1338 fr.block[b].sub[0].nr = nr[b];
1340 fr.block[b].sub[0].type = xdr_datatype_float;
1341 fr.block[b].sub[0].fval = block[b];
1343 fr.block[b].sub[0].type = xdr_datatype_double;
1344 fr.block[b].sub[0].dval = block[b];
1348 /* check for disre block & fill it. */
1353 add_blocks_enxframe(&fr, fr.nblock);
1355 add_subblocks_enxblock(&(fr.block[db]), 2);
1356 fr.block[db].id = enxDISRE;
1357 fr.block[db].sub[0].nr = ndisre;
1358 fr.block[db].sub[1].nr = ndisre;
1360 fr.block[db].sub[0].type = xdr_datatype_float;
1361 fr.block[db].sub[1].type = xdr_datatype_float;
1362 fr.block[db].sub[0].fval = disre_rt;
1363 fr.block[db].sub[1].fval = disre_rm3tav;
1365 fr.block[db].sub[0].type = xdr_datatype_double;
1366 fr.block[db].sub[1].type = xdr_datatype_double;
1367 fr.block[db].sub[0].dval = disre_rt;
1368 fr.block[db].sub[1].dval = disre_rm3tav;
1371 /* here we can put new-style blocks */
1373 /* Free energy perturbation blocks */
1376 mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1379 /* we can now free & reset the data in the blocks */
1382 mde_delta_h_coll_reset(md->dhc);
1385 /* AWH bias blocks. */
1386 if (awh != nullptr) // TODO: add boolean in t_mdebin. See in mdebin.h.
1388 awh->writeToEnergyFrame(step, &fr);
1391 /* do the actual I/O */
1392 do_enx(fp_ene, &fr);
1395 /* We have stored the sums, so reset the sum history */
1396 reset_ebin_sums(md->ebin);
1404 pprint(log, "A V E R A G E S", md);
1410 pprint(log, "R M S - F L U C T U A T I O N S", md);
1414 gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1419 for (i = 0; i < opts->ngtc; i++)
1421 if (opts->annealing[i] != eannNO)
1423 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1424 *(groups->grpname[groups->grps[egcTC].nm_ind[i]]),
1428 if (mode == eprNORMAL && fcd->orires.nr > 0)
1430 print_orires_log(log, &(fcd->orires));
1432 fprintf(log, " Energies (%s)\n", unit_energy);
1433 pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1436 if (mode == eprAVER)
1440 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1446 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1447 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1449 fprintf(log, " Force Virial (%s)\n", unit_energy);
1450 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1453 fprintf(log, " Total Virial (%s)\n", unit_energy);
1454 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1456 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1457 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1461 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1462 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1468 if (md->print_grpnms == nullptr)
1470 snew(md->print_grpnms, md->nE);
1472 for (i = 0; (i < md->nEg); i++)
1474 ni = groups->grps[egcENER].nm_ind[i];
1475 for (j = i; (j < md->nEg); j++)
1477 nj = groups->grps[egcENER].nm_ind[j];
1478 sprintf(buf, "%s-%s", *(groups->grpname[ni]),
1479 *(groups->grpname[nj]));
1480 md->print_grpnms[n++] = gmx_strdup(buf);
1484 sprintf(buf, "Epot (%s)", unit_energy);
1485 fprintf(log, "%15s ", buf);
1486 for (i = 0; (i < egNR); i++)
1490 fprintf(log, "%12s ", egrp_nm[i]);
1494 for (i = 0; (i < md->nE); i++)
1496 fprintf(log, "%15s", md->print_grpnms[i]);
1497 pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1504 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1509 fprintf(log, "%15s %12s %12s %12s\n",
1510 "Group", "Ux", "Uy", "Uz");
1511 for (i = 0; (i < md->nU); i++)
1513 ni = groups->grps[egcACC].nm_ind[i];
1514 fprintf(log, "%15s", *groups->grpname[ni]);
1515 pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1524 void update_energyhistory(energyhistory_t * enerhist, const t_mdebin * mdebin)
1526 const t_ebin * const ebin = mdebin->ebin;
1528 enerhist->nsteps = ebin->nsteps;
1529 enerhist->nsum = ebin->nsum;
1530 enerhist->nsteps_sim = ebin->nsteps_sim;
1531 enerhist->nsum_sim = ebin->nsum_sim;
1535 /* This will only actually resize the first time */
1536 enerhist->ener_ave.resize(ebin->nener);
1537 enerhist->ener_sum.resize(ebin->nener);
1539 for (int i = 0; i < ebin->nener; i++)
1541 enerhist->ener_ave[i] = ebin->e[i].eav;
1542 enerhist->ener_sum[i] = ebin->e[i].esum;
1546 if (ebin->nsum_sim > 0)
1548 /* This will only actually resize the first time */
1549 enerhist->ener_sum_sim.resize(ebin->nener);
1551 for (int i = 0; i < ebin->nener; i++)
1553 enerhist->ener_sum_sim[i] = ebin->e_sim[i].esum;
1558 mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1562 void restore_energyhistory_from_state(t_mdebin * mdebin,
1563 const energyhistory_t * enerhist)
1565 unsigned int nener = static_cast<unsigned int>(mdebin->ebin->nener);
1567 if ((enerhist->nsum > 0 && nener != enerhist->ener_sum.size()) ||
1568 (enerhist->nsum_sim > 0 && nener != enerhist->ener_sum_sim.size()))
1570 gmx_fatal(FARGS, "Mismatch between number of energies in run input (%d) and checkpoint file (%zu or %zu).",
1571 nener, enerhist->ener_sum.size(), enerhist->ener_sum_sim.size());
1574 mdebin->ebin->nsteps = enerhist->nsteps;
1575 mdebin->ebin->nsum = enerhist->nsum;
1576 mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1577 mdebin->ebin->nsum_sim = enerhist->nsum_sim;
1579 for (int i = 0; i < mdebin->ebin->nener; i++)
1581 mdebin->ebin->e[i].eav =
1582 (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1583 mdebin->ebin->e[i].esum =
1584 (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1585 mdebin->ebin->e_sim[i].esum =
1586 (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1590 mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist->deltaHForeignLambdas.get());