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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
57 #include "mtop_util.h"
61 #include "mdebin_bar.h"
64 static const char *conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
66 static const char *boxs_nm[] = { "Box-X", "Box-Y", "Box-Z" };
68 static const char *tricl_boxs_nm[] = {
69 "Box-XX", "Box-YY", "Box-ZZ",
70 "Box-YX", "Box-ZX", "Box-ZY"
73 static const char *vol_nm[] = { "Volume" };
75 static const char *dens_nm[] = {"Density" };
77 static const char *pv_nm[] = {"pV" };
79 static const char *enthalpy_nm[] = {"Enthalpy" };
81 static const char *boxvel_nm[] = {
82 "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
83 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
86 #define NBOXS asize(boxs_nm)
87 #define NTRICLBOXS asize(tricl_boxs_nm)
89 t_mdebin *init_mdebin(ener_file_t fp_ene,
90 const gmx_mtop_t *mtop,
94 const char *ener_nm[F_NRE];
95 static const char *vir_nm[] = {
96 "Vir-XX", "Vir-XY", "Vir-XZ",
97 "Vir-YX", "Vir-YY", "Vir-YZ",
98 "Vir-ZX", "Vir-ZY", "Vir-ZZ"
100 static const char *sv_nm[] = {
101 "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
102 "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
103 "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ"
105 static const char *fv_nm[] = {
106 "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
107 "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
108 "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ"
110 static const char *pres_nm[] = {
111 "Pres-XX", "Pres-XY", "Pres-XZ",
112 "Pres-YX", "Pres-YY", "Pres-YZ",
113 "Pres-ZX", "Pres-ZY", "Pres-ZZ"
115 static const char *surft_nm[] = {
118 static const char *mu_nm[] = {
119 "Mu-X", "Mu-Y", "Mu-Z"
121 static const char *vcos_nm[] = {
124 static const char *visc_nm[] = {
127 static const char *baro_nm[] = {
132 const gmx_groups_t *groups;
137 int i, j, ni, nj, n, nh, k, kk, ncon, nset;
138 gmx_bool bBHAM, bNoseHoover, b14;
147 if (EI_DYNAMICS(ir->eI))
149 md->delta_t = ir->delta_t;
156 groups = &mtop->groups;
158 bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
159 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
160 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
162 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
163 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
164 md->bConstr = (ncon > 0 || nset > 0);
165 md->bConstrVir = FALSE;
168 if (ncon > 0 && ir->eConstrAlg == econtLINCS)
179 md->bConstrVir = (getenv("GMX_CONSTRAINTVIR") != NULL);
186 /* Energy monitoring */
187 for (i = 0; i < egNR; i++)
189 md->bEInd[i] = FALSE;
192 /* Even though the OpenMM build has moved to contrib, it's not
193 * practical to move/remove this code fragment, because of the
194 * fundamental mess that is the GROMACS library structure. */
196 for (i = 0; i < F_NRE; i++)
198 md->bEner[i] = FALSE;
201 md->bEner[i] = !bBHAM;
203 else if (i == F_BHAM)
205 md->bEner[i] = bBHAM;
209 md->bEner[i] = ir->bQMMM;
211 else if (i == F_COUL_LR)
213 md->bEner[i] = (ir->rcoulomb > ir->rlist);
215 else if (i == F_LJ_LR)
217 md->bEner[i] = (!bBHAM && ir->rvdw > ir->rlist);
219 else if (i == F_BHAM_LR)
221 md->bEner[i] = (bBHAM && ir->rvdw > ir->rlist);
223 else if (i == F_RF_EXCL)
225 md->bEner[i] = (EEL_RF(ir->coulombtype) && ir->coulombtype != eelRF_NEC && ir->cutoff_scheme == ecutsGROUP);
227 else if (i == F_COUL_RECIP)
229 md->bEner[i] = EEL_FULL(ir->coulombtype);
231 else if (i == F_LJ14)
235 else if (i == F_COUL14)
239 else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
241 md->bEner[i] = FALSE;
243 else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
244 (i == F_DVDL_VDW && ir->fepvals->separate_dvdl[efptVDW]) ||
245 (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
246 (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
247 (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
248 (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
250 md->bEner[i] = (ir->efep != efepNO);
252 else if ((interaction_function[i].flags & IF_VSITE) ||
253 (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
255 md->bEner[i] = FALSE;
257 else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES) || (i == F_EQM))
261 else if ((i == F_GBPOL) && ir->implicit_solvent == eisGBSA)
265 else if ((i == F_NPSOLVATION) && ir->implicit_solvent == eisGBSA && (ir->sa_algorithm != esaNO))
269 else if ((i == F_GB12) || (i == F_GB13) || (i == F_GB14))
271 md->bEner[i] = FALSE;
273 else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
275 md->bEner[i] = EI_DYNAMICS(ir->eI);
277 else if (i == F_DISPCORR || i == F_PDISPCORR)
279 md->bEner[i] = (ir->eDispCorr != edispcNO);
281 else if (i == F_DISRESVIOL)
283 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
285 else if (i == F_ORIRESDEV)
287 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
289 else if (i == F_CONNBONDS)
291 md->bEner[i] = FALSE;
293 else if (i == F_COM_PULL)
295 md->bEner[i] = (ir->ePull == epullUMBRELLA || ir->ePull == epullCONST_F || ir->bRot);
297 else if (i == F_ECONSERVED)
299 md->bEner[i] = ((ir->etc == etcNOSEHOOVER || ir->etc == etcVRESCALE) &&
300 (ir->epc == epcNO || ir->epc == epcMTTK));
304 md->bEner[i] = (gmx_mtop_ftype_count(mtop, i) > 0);
308 /* OpenMM always produces only the following 4 energy terms */
309 md->bEner[F_EPOT] = TRUE;
310 md->bEner[F_EKIN] = TRUE;
311 md->bEner[F_ETOT] = TRUE;
312 md->bEner[F_TEMP] = TRUE;
315 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
316 if (ir->bAdress && !debug)
318 for (i = 0; i < F_NRE; i++)
320 md->bEner[i] = FALSE;
337 for (i = 0; i < F_NRE; i++)
341 ener_nm[md->f_nre] = interaction_function[i].longname;
347 md->bDiagPres = !TRICLINIC(ir->ref_p);
348 md->ref_p = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
349 md->bTricl = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
350 md->bDynBox = DYNAMIC_BOX(*ir);
352 md->bNHC_trotter = IR_NVT_TROTTER(ir);
353 md->bPrintNHChains = ir->bPrintNHChains;
354 md->bMTTK = (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir));
355 md->bMu = NEED_MUTOT(*ir);
357 md->ebin = mk_ebin();
358 /* Pass NULL for unit to let get_ebin_space determine the units
359 * for interaction_function[i].longname
361 md->ie = get_ebin_space(md->ebin, md->f_nre, ener_nm, NULL);
364 /* This should be called directly after the call for md->ie,
365 * such that md->iconrmsd follows directly in the list.
367 md->iconrmsd = get_ebin_space(md->ebin, md->nCrmsd, conrmsd_nm, "");
371 md->ib = get_ebin_space(md->ebin,
372 md->bTricl ? NTRICLBOXS : NBOXS,
373 md->bTricl ? tricl_boxs_nm : boxs_nm,
375 md->ivol = get_ebin_space(md->ebin, 1, vol_nm, unit_volume);
376 md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
379 md->ipv = get_ebin_space(md->ebin, 1, pv_nm, unit_energy);
380 md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm, unit_energy);
385 md->isvir = get_ebin_space(md->ebin, asize(sv_nm), sv_nm, unit_energy);
386 md->ifvir = get_ebin_space(md->ebin, asize(fv_nm), fv_nm, unit_energy);
390 md->ivir = get_ebin_space(md->ebin, asize(vir_nm), vir_nm, unit_energy);
394 md->ipres = get_ebin_space(md->ebin, asize(pres_nm), pres_nm, unit_pres_bar);
398 md->isurft = get_ebin_space(md->ebin, asize(surft_nm), surft_nm,
401 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
403 md->ipc = get_ebin_space(md->ebin, md->bTricl ? 6 : 3,
404 boxvel_nm, unit_vel);
408 md->imu = get_ebin_space(md->ebin, asize(mu_nm), mu_nm, unit_dipole_D);
410 if (ir->cos_accel != 0)
412 md->ivcos = get_ebin_space(md->ebin, asize(vcos_nm), vcos_nm, unit_vel);
413 md->ivisc = get_ebin_space(md->ebin, asize(visc_nm), visc_nm,
417 /* Energy monitoring */
418 for (i = 0; i < egNR; i++)
420 md->bEInd[i] = FALSE;
422 md->bEInd[egCOULSR] = TRUE;
423 md->bEInd[egLJSR ] = TRUE;
425 if (ir->rcoulomb > ir->rlist)
427 md->bEInd[egCOULLR] = TRUE;
431 if (ir->rvdw > ir->rlist)
433 md->bEInd[egLJLR] = TRUE;
438 md->bEInd[egLJSR] = FALSE;
439 md->bEInd[egBHAMSR] = TRUE;
440 if (ir->rvdw > ir->rlist)
442 md->bEInd[egBHAMLR] = TRUE;
447 md->bEInd[egLJ14] = TRUE;
448 md->bEInd[egCOUL14] = TRUE;
451 for (i = 0; (i < egNR); i++)
459 n = groups->grps[egcENER].nr;
460 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
463 /*standard simulation*/
465 md->nE = (n*(n+1))/2;
469 /*AdResS simulation*/
476 snew(md->igrp, md->nE);
481 for (k = 0; (k < md->nEc); k++)
483 snew(gnm[k], STRLEN);
485 for (i = 0; (i < groups->grps[egcENER].nr); i++)
487 ni = groups->grps[egcENER].nm_ind[i];
488 for (j = i; (j < groups->grps[egcENER].nr); j++)
490 nj = groups->grps[egcENER].nm_ind[j];
491 for (k = kk = 0; (k < egNR); k++)
495 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
496 *(groups->grpname[ni]), *(groups->grpname[nj]));
500 md->igrp[n] = get_ebin_space(md->ebin, md->nEc,
501 (const char **)gnm, unit_energy);
505 for (k = 0; (k < md->nEc); k++)
513 gmx_incons("Number of energy terms wrong");
517 md->nTC = groups->grps[egcTC].nr;
518 md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
521 md->nTCP = 1; /* assume only one possible coupling system for barostat
528 if (md->etc == etcNOSEHOOVER)
530 if (md->bNHC_trotter)
532 md->mde_n = 2*md->nNHC*md->nTC;
536 md->mde_n = 2*md->nTC;
538 if (md->epc == epcMTTK)
540 md->mdeb_n = 2*md->nNHC*md->nTCP;
549 snew(md->tmp_r, md->mde_n);
550 snew(md->tmp_v, md->mde_n);
551 snew(md->grpnms, md->mde_n);
554 for (i = 0; (i < md->nTC); i++)
556 ni = groups->grps[egcTC].nm_ind[i];
557 sprintf(buf, "T-%s", *(groups->grpname[ni]));
558 grpnms[i] = strdup(buf);
560 md->itemp = get_ebin_space(md->ebin, md->nTC, (const char **)grpnms,
563 if (md->etc == etcNOSEHOOVER)
565 if (md->bPrintNHChains)
567 if (md->bNHC_trotter)
569 for (i = 0; (i < md->nTC); i++)
571 ni = groups->grps[egcTC].nm_ind[i];
572 bufi = *(groups->grpname[ni]);
573 for (j = 0; (j < md->nNHC); j++)
575 sprintf(buf, "Xi-%d-%s", j, bufi);
576 grpnms[2*(i*md->nNHC+j)] = strdup(buf);
577 sprintf(buf, "vXi-%d-%s", j, bufi);
578 grpnms[2*(i*md->nNHC+j)+1] = strdup(buf);
581 md->itc = get_ebin_space(md->ebin, md->mde_n,
582 (const char **)grpnms, unit_invtime);
585 for (i = 0; (i < md->nTCP); i++)
587 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
588 for (j = 0; (j < md->nNHC); j++)
590 sprintf(buf, "Xi-%d-%s", j, bufi);
591 grpnms[2*(i*md->nNHC+j)] = strdup(buf);
592 sprintf(buf, "vXi-%d-%s", j, bufi);
593 grpnms[2*(i*md->nNHC+j)+1] = strdup(buf);
596 md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
597 (const char **)grpnms, unit_invtime);
602 for (i = 0; (i < md->nTC); i++)
604 ni = groups->grps[egcTC].nm_ind[i];
605 bufi = *(groups->grpname[ni]);
606 sprintf(buf, "Xi-%s", bufi);
607 grpnms[2*i] = strdup(buf);
608 sprintf(buf, "vXi-%s", bufi);
609 grpnms[2*i+1] = strdup(buf);
611 md->itc = get_ebin_space(md->ebin, md->mde_n,
612 (const char **)grpnms, unit_invtime);
616 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
617 md->etc == etcVRESCALE)
619 for (i = 0; (i < md->nTC); i++)
621 ni = groups->grps[egcTC].nm_ind[i];
622 sprintf(buf, "Lamb-%s", *(groups->grpname[ni]));
623 grpnms[i] = strdup(buf);
625 md->itc = get_ebin_space(md->ebin, md->mde_n, (const char **)grpnms, "");
631 md->nU = groups->grps[egcACC].nr;
634 snew(grpnms, 3*md->nU);
635 for (i = 0; (i < md->nU); i++)
637 ni = groups->grps[egcACC].nm_ind[i];
638 sprintf(buf, "Ux-%s", *(groups->grpname[ni]));
639 grpnms[3*i+XX] = strdup(buf);
640 sprintf(buf, "Uy-%s", *(groups->grpname[ni]));
641 grpnms[3*i+YY] = strdup(buf);
642 sprintf(buf, "Uz-%s", *(groups->grpname[ni]));
643 grpnms[3*i+ZZ] = strdup(buf);
645 md->iu = get_ebin_space(md->ebin, 3*md->nU, (const char **)grpnms, unit_vel);
651 do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
654 md->print_grpnms = NULL;
656 /* check whether we're going to write dh histograms */
658 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
660 /* Currently dh histograms are only written with dynamics */
661 if (EI_DYNAMICS(ir->eI))
665 mde_delta_h_coll_init(md->dhc, ir);
668 snew(md->dE, ir->fepvals->n_lambda);
672 md->fp_dhdl = fp_dhdl;
673 snew(md->dE, ir->fepvals->n_lambda);
678 snew(md->temperatures, ir->fepvals->n_lambda);
679 for (i = 0; i < ir->fepvals->n_lambda; i++)
681 md->temperatures[i] = ir->simtempvals->temperatures[i];
687 /* print a lambda vector to a string
688 fep = the inputrec's FEP input data
689 i = the index of the lambda vector
690 get_native_lambda = whether to print the native lambda
691 get_names = whether to print the names rather than the values
692 str = the pre-allocated string buffer to print to. */
693 static void print_lambda_vector(t_lambda *fep, int i,
694 gmx_bool get_native_lambda, gmx_bool get_names,
701 for (j = 0; j < efptNR; j++)
703 if (fep->separate_dvdl[j])
708 str[0] = 0; /* reset the string */
711 str += sprintf(str, "("); /* set the opening parenthesis*/
713 for (j = 0; j < efptNR; j++)
715 if (fep->separate_dvdl[j])
720 if (get_native_lambda && fep->init_lambda >= 0)
722 str += sprintf(str, "%.4f", fep->init_lambda);
726 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
731 str += sprintf(str, "%s", efpt_singular_names[j]);
733 /* print comma for the next item */
736 str += sprintf(str, ", ");
743 /* and add the closing parenthesis */
744 str += sprintf(str, ")");
749 extern FILE *open_dhdl(const char *filename, const t_inputrec *ir,
750 const output_env_t oenv)
753 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
754 *lambdastate = "\\lambda state", *remain = "remaining";
755 char title[STRLEN], label_x[STRLEN], label_y[STRLEN];
756 int i, np, nps, nsets, nsets_de, nsetsbegin;
757 int n_lambda_terms = 0;
758 t_lambda *fep = ir->fepvals; /* for simplicity */
759 t_expanded *expand = ir->expandedvals;
761 char buf[STRLEN], lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
767 gmx_bool write_pV = FALSE;
769 /* count the number of different lambda terms */
770 for (i = 0; i < efptNR; i++)
772 if (fep->separate_dvdl[i])
778 if (fep->n_lambda == 0)
780 sprintf(title, "%s", dhdl);
781 sprintf(label_x, "Time (ps)");
782 sprintf(label_y, "%s (%s %s)",
783 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
787 sprintf(title, "%s and %s", dhdl, deltag);
788 sprintf(label_x, "Time (ps)");
789 sprintf(label_y, "%s and %s (%s %s)",
790 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
792 fp = gmx_fio_fopen(filename, "w+");
793 xvgr_header(fp, title, label_x, label_y, exvggtXNY, oenv);
797 bufplace = sprintf(buf, "T = %g (K) ",
800 if (ir->efep != efepSLOWGROWTH)
802 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
804 /* compatibility output */
805 sprintf(&(buf[bufplace]), "%s = %.4f", lambda, fep->init_lambda);
809 print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
811 print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
813 sprintf(&(buf[bufplace]), "%s %d: %s = %s",
814 lambdastate, fep->init_fep_state,
815 lambda_name_str, lambda_vec_str);
818 xvgr_subtitle(fp, buf, oenv);
822 if (fep->dhdl_derivatives == edhdlderivativesYES)
824 nsets_dhdl = n_lambda_terms;
826 /* count the number of delta_g states */
827 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
829 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
831 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
833 nsets += 1; /*add fep state for expanded ensemble */
836 if (fep->bPrintEnergy)
838 nsets += 1; /* add energy to the dhdl as well */
842 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
844 nsetsextend += 1; /* for PV term, other terms possible if required for
845 the reduced potential (only needed with foreign
846 lambda, and only output when init_lambda is not
847 set in order to maintain compatibility of the
851 snew(setname, nsetsextend);
853 if (expand->elmcmove > elmcmoveNO)
855 /* state for the fep_vals, if we have alchemical sampling */
856 sprintf(buf, "%s", "Thermodynamic state");
857 setname[s] = strdup(buf);
861 if (fep->bPrintEnergy)
863 sprintf(buf, "%s (%s)", "Energy", unit_energy);
864 setname[s] = strdup(buf);
868 if (fep->dhdl_derivatives == edhdlderivativesYES)
870 for (i = 0; i < efptNR; i++)
872 if (fep->separate_dvdl[i])
875 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
877 /* compatibility output */
878 sprintf(buf, "%s %s %.4f", dhdl, lambda, fep->init_lambda);
882 double lam = fep->init_lambda;
883 if (fep->init_lambda < 0)
885 lam = fep->all_lambda[i][fep->init_fep_state];
887 sprintf(buf, "%s %s = %.4f", dhdl, efpt_singular_names[i],
890 setname[s] = strdup(buf);
896 if (fep->n_lambda > 0)
898 /* g_bar has to determine the lambda values used in this simulation
899 * from this xvg legend.
902 if (expand->elmcmove > elmcmoveNO)
904 nsetsbegin = 1; /* for including the expanded ensemble */
911 if (fep->bPrintEnergy)
915 nsetsbegin += nsets_dhdl;
917 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
919 print_lambda_vector(fep, i, FALSE, FALSE, lambda_vec_str);
920 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
922 /* for compatible dhdl.xvg files */
923 nps = sprintf(buf, "%s %s %s", deltag, lambda, lambda_vec_str);
927 nps = sprintf(buf, "%s %s to %s", deltag, lambda, lambda_vec_str);
932 /* print the temperature for this state if doing simulated annealing */
933 sprintf(&buf[nps], "T = %g (%s)",
934 ir->simtempvals->temperatures[s-(nsetsbegin)],
937 setname[s] = strdup(buf);
942 np = sprintf(buf, "pV (%s)", unit_energy);
943 setname[nsetsextend-1] = strdup(buf); /* the first entry after
947 xvgr_legend(fp, nsetsextend, (const char **)setname, oenv);
949 for (s = 0; s < nsetsextend; s++)
959 static void copy_energy(t_mdebin *md, real e[], real ecpy[])
963 for (i = j = 0; (i < F_NRE); i++)
972 gmx_incons("Number of energy terms wrong");
976 void upd_mdebin(t_mdebin *md,
981 gmx_enerdata_t *enerd,
990 gmx_ekindata_t *ekind,
994 int i, j, k, kk, m, n, gid;
995 real crmsd[2], tmp6[6];
996 real bs[NTRICLBOXS], vol, dens, pv, enthalpy;
999 double store_dhdl[efptNR];
1000 real store_energy = 0;
1003 /* Do NOT use the box in the state variable, but the separate box provided
1004 * as an argument. This is because we sometimes need to write the box from
1005 * the last timestep to match the trajectory frames.
1007 copy_energy(md, enerd->term, ecopy);
1008 add_ebin(md->ebin, md->ie, md->f_nre, ecopy, bSum);
1011 crmsd[0] = constr_rmsd(constr, FALSE);
1014 crmsd[1] = constr_rmsd(constr, TRUE);
1016 add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
1023 bs[0] = box[XX][XX];
1024 bs[1] = box[YY][YY];
1025 bs[2] = box[ZZ][ZZ];
1026 bs[3] = box[YY][XX];
1027 bs[4] = box[ZZ][XX];
1028 bs[5] = box[ZZ][YY];
1033 bs[0] = box[XX][XX];
1034 bs[1] = box[YY][YY];
1035 bs[2] = box[ZZ][ZZ];
1038 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1039 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1040 add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
1041 add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
1042 add_ebin(md->ebin, md->idens, 1, &dens, bSum);
1046 /* This is pV (in kJ/mol). The pressure is the reference pressure,
1047 not the instantaneous pressure */
1048 pv = vol*md->ref_p/PRESFAC;
1050 add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
1051 enthalpy = pv + enerd->term[F_ETOT];
1052 add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
1057 add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
1058 add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
1062 add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
1066 add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
1070 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
1071 add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
1073 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
1075 tmp6[0] = state->boxv[XX][XX];
1076 tmp6[1] = state->boxv[YY][YY];
1077 tmp6[2] = state->boxv[ZZ][ZZ];
1078 tmp6[3] = state->boxv[YY][XX];
1079 tmp6[4] = state->boxv[ZZ][XX];
1080 tmp6[5] = state->boxv[ZZ][YY];
1081 add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1085 add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1087 if (ekind && ekind->cosacc.cos_accel != 0)
1089 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1090 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1091 add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1092 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1093 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1094 *dens*vol*sqr(box[ZZ][ZZ]*NANO/(2*M_PI)));
1095 add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1100 for (i = 0; (i < md->nEg); i++)
1102 for (j = i; (j < md->nEg); j++)
1104 gid = GID(i, j, md->nEg);
1105 for (k = kk = 0; (k < egNR); k++)
1109 eee[kk++] = enerd->grpp.ener[k][gid];
1112 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1120 for (i = 0; (i < md->nTC); i++)
1122 md->tmp_r[i] = ekind->tcstat[i].T;
1124 add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1126 if (md->etc == etcNOSEHOOVER)
1128 /* whether to print Nose-Hoover chains: */
1129 if (md->bPrintNHChains)
1131 if (md->bNHC_trotter)
1133 for (i = 0; (i < md->nTC); i++)
1135 for (j = 0; j < md->nNHC; j++)
1138 md->tmp_r[2*k] = state->nosehoover_xi[k];
1139 md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1142 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1146 for (i = 0; (i < md->nTCP); i++)
1148 for (j = 0; j < md->nNHC; j++)
1151 md->tmp_r[2*k] = state->nhpres_xi[k];
1152 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1155 add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1160 for (i = 0; (i < md->nTC); i++)
1162 md->tmp_r[2*i] = state->nosehoover_xi[i];
1163 md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1165 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1169 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1170 md->etc == etcVRESCALE)
1172 for (i = 0; (i < md->nTC); i++)
1174 md->tmp_r[i] = ekind->tcstat[i].lambda;
1176 add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1180 if (ekind && md->nU > 1)
1182 for (i = 0; (i < md->nU); i++)
1184 copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1186 add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1189 ebin_increase_count(md->ebin, bSum);
1191 /* BAR + thermodynamic integration values */
1192 if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1194 for (i = 0; i < enerd->n_lambda-1; i++)
1196 /* zero for simulated tempering */
1197 md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1198 if (md->temperatures != NULL)
1200 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1201 /* is this even useful to have at all? */
1202 md->dE[i] += (md->temperatures[i]/
1203 md->temperatures[state->fep_state]-1.0)*
1204 enerd->term[F_EKIN];
1210 fprintf(md->fp_dhdl, "%.4f", time);
1211 /* the current free energy state */
1213 /* print the current state if we are doing expanded ensemble */
1214 if (expand->elmcmove > elmcmoveNO)
1216 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1218 /* total energy (for if the temperature changes */
1219 if (fep->bPrintEnergy)
1221 store_energy = enerd->term[F_ETOT];
1222 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1225 if (fep->dhdl_derivatives == edhdlderivativesYES)
1227 for (i = 0; i < efptNR; i++)
1229 if (fep->separate_dvdl[i])
1231 /* assumes F_DVDL is first */
1232 fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1236 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1238 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1240 if ((md->epc != epcNO) &&
1241 (enerd->n_lambda > 0) &&
1242 (fep->init_lambda < 0))
1244 fprintf(md->fp_dhdl, " %#.8g", pv); /* PV term only needed when
1245 there are alternate state
1246 lambda and we're not in
1247 compatibility mode */
1249 fprintf(md->fp_dhdl, "\n");
1250 /* and the binary free energy output */
1252 if (md->dhc && bDoDHDL)
1255 for (i = 0; i < efptNR; i++)
1257 if (fep->separate_dvdl[i])
1259 /* assumes F_DVDL is first */
1260 store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1264 store_energy = enerd->term[F_ETOT];
1265 /* store_dh is dE */
1266 mde_delta_h_coll_add_dh(md->dhc,
1267 (double)state->fep_state,
1271 md->dE + fep->lambda_start_n,
1278 void upd_mdebin_step(t_mdebin *md)
1280 ebin_increase_count(md->ebin, FALSE);
1283 static void npr(FILE *log, int n, char c)
1285 for (; (n > 0); n--)
1287 fprintf(log, "%c", c);
1291 static void pprint(FILE *log, const char *s, t_mdebin *md)
1295 char buf1[22], buf2[22];
1298 fprintf(log, "\t<====== ");
1299 npr(log, slen, CHAR);
1300 fprintf(log, " ==>\n");
1301 fprintf(log, "\t<==== %s ====>\n", s);
1302 fprintf(log, "\t<== ");
1303 npr(log, slen, CHAR);
1304 fprintf(log, " ======>\n\n");
1306 fprintf(log, "\tStatistics over %s steps using %s frames\n",
1307 gmx_step_str(md->ebin->nsteps_sim, buf1),
1308 gmx_step_str(md->ebin->nsum_sim, buf2));
1312 void print_ebin_header(FILE *log, gmx_large_int_t steps, double time, real lambda)
1316 fprintf(log, " %12s %12s %12s\n"
1317 " %12s %12.5f %12.5f\n\n",
1318 "Step", "Time", "Lambda", gmx_step_str(steps, buf), time, lambda);
1321 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1323 gmx_large_int_t step, double time,
1324 int mode, gmx_bool bCompact,
1325 t_mdebin *md, t_fcdata *fcd,
1326 gmx_groups_t *groups, t_grpopts *opts)
1328 /*static char **grpnms=NULL;*/
1330 int i, j, n, ni, nj, ndr, nor, b;
1332 real *disre_rm3tav, *disre_rt;
1334 /* these are for the old-style blocks (1 subblock, only reals), because
1335 there can be only one per ID for these */
1340 /* temporary arrays for the lambda values to write out */
1341 double enxlambda_data[2];
1351 fr.nsteps = md->ebin->nsteps;
1352 fr.dt = md->delta_t;
1353 fr.nsum = md->ebin->nsum;
1354 fr.nre = (bEne) ? md->ebin->nener : 0;
1355 fr.ener = md->ebin->e;
1356 ndisre = bDR ? fcd->disres.npair : 0;
1357 disre_rm3tav = fcd->disres.rm3tav;
1358 disre_rt = fcd->disres.rt;
1359 /* Optional additional old-style (real-only) blocks. */
1360 for (i = 0; i < enxNR; i++)
1364 if (fcd->orires.nr > 0 && bOR)
1366 diagonalize_orires_tensors(&(fcd->orires));
1367 nr[enxOR] = fcd->orires.nr;
1368 block[enxOR] = fcd->orires.otav;
1370 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1372 block[enxORI] = fcd->orires.oinsl;
1373 id[enxORI] = enxORI;
1374 nr[enxORT] = fcd->orires.nex*12;
1375 block[enxORT] = fcd->orires.eig;
1376 id[enxORT] = enxORT;
1379 /* whether we are going to wrte anything out: */
1380 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1383 /* the old-style blocks go first */
1385 for (i = 0; i < enxNR; i++)
1392 add_blocks_enxframe(&fr, fr.nblock);
1393 for (b = 0; b < fr.nblock; b++)
1395 add_subblocks_enxblock(&(fr.block[b]), 1);
1396 fr.block[b].id = id[b];
1397 fr.block[b].sub[0].nr = nr[b];
1399 fr.block[b].sub[0].type = xdr_datatype_float;
1400 fr.block[b].sub[0].fval = block[b];
1402 fr.block[b].sub[0].type = xdr_datatype_double;
1403 fr.block[b].sub[0].dval = block[b];
1407 /* check for disre block & fill it. */
1412 add_blocks_enxframe(&fr, fr.nblock);
1414 add_subblocks_enxblock(&(fr.block[db]), 2);
1415 fr.block[db].id = enxDISRE;
1416 fr.block[db].sub[0].nr = ndisre;
1417 fr.block[db].sub[1].nr = ndisre;
1419 fr.block[db].sub[0].type = xdr_datatype_float;
1420 fr.block[db].sub[1].type = xdr_datatype_float;
1421 fr.block[db].sub[0].fval = disre_rt;
1422 fr.block[db].sub[1].fval = disre_rm3tav;
1424 fr.block[db].sub[0].type = xdr_datatype_double;
1425 fr.block[db].sub[1].type = xdr_datatype_double;
1426 fr.block[db].sub[0].dval = disre_rt;
1427 fr.block[db].sub[1].dval = disre_rm3tav;
1430 /* here we can put new-style blocks */
1432 /* Free energy perturbation blocks */
1435 mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1438 /* we can now free & reset the data in the blocks */
1441 mde_delta_h_coll_reset(md->dhc);
1444 /* do the actual I/O */
1445 do_enx(fp_ene, &fr);
1446 gmx_fio_check_file_position(enx_file_pointer(fp_ene));
1449 /* We have stored the sums, so reset the sum history */
1450 reset_ebin_sums(md->ebin);
1458 pprint(log, "A V E R A G E S", md);
1464 pprint(log, "R M S - F L U C T U A T I O N S", md);
1468 gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1473 for (i = 0; i < opts->ngtc; i++)
1475 if (opts->annealing[i] != eannNO)
1477 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1478 *(groups->grpname[groups->grps[egcTC].nm_ind[i]]),
1482 if (mode == eprNORMAL && fcd->orires.nr > 0)
1484 print_orires_log(log, &(fcd->orires));
1486 fprintf(log, " Energies (%s)\n", unit_energy);
1487 pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1494 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1500 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1501 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1503 fprintf(log, " Force Virial (%s)\n", unit_energy);
1504 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1509 fprintf(log, " Total Virial (%s)\n", unit_energy);
1510 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1515 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1516 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1521 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1522 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1528 if (md->print_grpnms == NULL)
1530 snew(md->print_grpnms, md->nE);
1532 for (i = 0; (i < md->nEg); i++)
1534 ni = groups->grps[egcENER].nm_ind[i];
1535 for (j = i; (j < md->nEg); j++)
1537 nj = groups->grps[egcENER].nm_ind[j];
1538 sprintf(buf, "%s-%s", *(groups->grpname[ni]),
1539 *(groups->grpname[nj]));
1540 md->print_grpnms[n++] = strdup(buf);
1544 sprintf(buf, "Epot (%s)", unit_energy);
1545 fprintf(log, "%15s ", buf);
1546 for (i = 0; (i < egNR); i++)
1550 fprintf(log, "%12s ", egrp_nm[i]);
1554 for (i = 0; (i < md->nE); i++)
1556 fprintf(log, "%15s", md->print_grpnms[i]);
1557 pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1564 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1569 fprintf(log, "%15s %12s %12s %12s\n",
1570 "Group", "Ux", "Uy", "Uz");
1571 for (i = 0; (i < md->nU); i++)
1573 ni = groups->grps[egcACC].nm_ind[i];
1574 fprintf(log, "%15s", *groups->grpname[ni]);
1575 pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1584 void update_energyhistory(energyhistory_t * enerhist, t_mdebin * mdebin)
1588 enerhist->nsteps = mdebin->ebin->nsteps;
1589 enerhist->nsum = mdebin->ebin->nsum;
1590 enerhist->nsteps_sim = mdebin->ebin->nsteps_sim;
1591 enerhist->nsum_sim = mdebin->ebin->nsum_sim;
1592 enerhist->nener = mdebin->ebin->nener;
1594 if (mdebin->ebin->nsum > 0)
1596 /* Check if we need to allocate first */
1597 if (enerhist->ener_ave == NULL)
1599 snew(enerhist->ener_ave, enerhist->nener);
1600 snew(enerhist->ener_sum, enerhist->nener);
1603 for (i = 0; i < enerhist->nener; i++)
1605 enerhist->ener_ave[i] = mdebin->ebin->e[i].eav;
1606 enerhist->ener_sum[i] = mdebin->ebin->e[i].esum;
1610 if (mdebin->ebin->nsum_sim > 0)
1612 /* Check if we need to allocate first */
1613 if (enerhist->ener_sum_sim == NULL)
1615 snew(enerhist->ener_sum_sim, enerhist->nener);
1618 for (i = 0; i < enerhist->nener; i++)
1620 enerhist->ener_sum_sim[i] = mdebin->ebin->e_sim[i].esum;
1625 mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1629 void restore_energyhistory_from_state(t_mdebin * mdebin,
1630 energyhistory_t * enerhist)
1634 if ((enerhist->nsum > 0 || enerhist->nsum_sim > 0) &&
1635 mdebin->ebin->nener != enerhist->nener)
1637 gmx_fatal(FARGS, "Mismatch between number of energies in run input (%d) and checkpoint file (%d).",
1638 mdebin->ebin->nener, enerhist->nener);
1641 mdebin->ebin->nsteps = enerhist->nsteps;
1642 mdebin->ebin->nsum = enerhist->nsum;
1643 mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1644 mdebin->ebin->nsum_sim = enerhist->nsum_sim;
1646 for (i = 0; i < mdebin->ebin->nener; i++)
1648 mdebin->ebin->e[i].eav =
1649 (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1650 mdebin->ebin->e[i].esum =
1651 (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1652 mdebin->ebin->e_sim[i].esum =
1653 (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1657 mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist);