1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
55 #include "mtop_util.h"
60 #include "mdebin_bar.h"
63 static const char *conrmsd_nm[] = { "Constr. rmsd", "Constr.2 rmsd" };
65 static const char *boxs_nm[] = { "Box-X", "Box-Y", "Box-Z" };
67 static const char *tricl_boxs_nm[] = {
68 "Box-XX", "Box-YY", "Box-ZZ",
69 "Box-YX", "Box-ZX", "Box-ZY"
72 static const char *vol_nm[] = { "Volume" };
74 static const char *dens_nm[] = {"Density" };
76 static const char *pv_nm[] = {"pV" };
78 static const char *enthalpy_nm[] = {"Enthalpy" };
80 static const char *boxvel_nm[] = {
81 "Box-Vel-XX", "Box-Vel-YY", "Box-Vel-ZZ",
82 "Box-Vel-YX", "Box-Vel-ZX", "Box-Vel-ZY"
85 #define NBOXS asize(boxs_nm)
86 #define NTRICLBOXS asize(tricl_boxs_nm)
88 t_mdebin *init_mdebin(ener_file_t fp_ene,
89 const gmx_mtop_t *mtop,
93 const char *ener_nm[F_NRE];
94 static const char *vir_nm[] = {
95 "Vir-XX", "Vir-XY", "Vir-XZ",
96 "Vir-YX", "Vir-YY", "Vir-YZ",
97 "Vir-ZX", "Vir-ZY", "Vir-ZZ"
99 static const char *sv_nm[] = {
100 "ShakeVir-XX", "ShakeVir-XY", "ShakeVir-XZ",
101 "ShakeVir-YX", "ShakeVir-YY", "ShakeVir-YZ",
102 "ShakeVir-ZX", "ShakeVir-ZY", "ShakeVir-ZZ"
104 static const char *fv_nm[] = {
105 "ForceVir-XX", "ForceVir-XY", "ForceVir-XZ",
106 "ForceVir-YX", "ForceVir-YY", "ForceVir-YZ",
107 "ForceVir-ZX", "ForceVir-ZY", "ForceVir-ZZ"
109 static const char *pres_nm[] = {
110 "Pres-XX", "Pres-XY", "Pres-XZ",
111 "Pres-YX", "Pres-YY", "Pres-YZ",
112 "Pres-ZX", "Pres-ZY", "Pres-ZZ"
114 static const char *surft_nm[] = {
117 static const char *mu_nm[] = {
118 "Mu-X", "Mu-Y", "Mu-Z"
120 static const char *vcos_nm[] = {
123 static const char *visc_nm[] = {
126 static const char *baro_nm[] = {
131 const gmx_groups_t *groups;
136 int i, j, ni, nj, n, nh, k, kk, ncon, nset;
137 gmx_bool bBHAM, bNoseHoover, b14;
146 if (EI_DYNAMICS(ir->eI))
148 md->delta_t = ir->delta_t;
155 groups = &mtop->groups;
157 bBHAM = (mtop->ffparams.functype[0] == F_BHAM);
158 b14 = (gmx_mtop_ftype_count(mtop, F_LJ14) > 0 ||
159 gmx_mtop_ftype_count(mtop, F_LJC14_Q) > 0);
161 ncon = gmx_mtop_ftype_count(mtop, F_CONSTR);
162 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
163 md->bConstr = (ncon > 0 || nset > 0);
164 md->bConstrVir = FALSE;
167 if (ncon > 0 && ir->eConstrAlg == econtLINCS)
178 md->bConstrVir = (getenv("GMX_CONSTRAINTVIR") != NULL);
185 /* Energy monitoring */
186 for (i = 0; i < egNR; i++)
188 md->bEInd[i] = FALSE;
191 /* Even though the OpenMM build has moved to contrib, it's not
192 * practical to move/remove this code fragment, because of the
193 * fundamental mess that is the GROMACS library structure. */
195 for (i = 0; i < F_NRE; i++)
197 md->bEner[i] = FALSE;
200 md->bEner[i] = !bBHAM;
202 else if (i == F_BHAM)
204 md->bEner[i] = bBHAM;
208 md->bEner[i] = ir->bQMMM;
210 else if (i == F_COUL_LR)
212 md->bEner[i] = (ir->rcoulomb > ir->rlist);
214 else if (i == F_LJ_LR)
216 md->bEner[i] = (!bBHAM && ir->rvdw > ir->rlist);
218 else if (i == F_BHAM_LR)
220 md->bEner[i] = (bBHAM && ir->rvdw > ir->rlist);
222 else if (i == F_RF_EXCL)
224 md->bEner[i] = (EEL_RF(ir->coulombtype) && ir->coulombtype != eelRF_NEC && ir->cutoff_scheme == ecutsGROUP);
226 else if (i == F_COUL_RECIP)
228 md->bEner[i] = EEL_FULL(ir->coulombtype);
230 else if (i == F_LJ14)
234 else if (i == F_COUL14)
238 else if (i == F_LJC14_Q || i == F_LJC_PAIRS_NB)
240 md->bEner[i] = FALSE;
242 else if ((i == F_DVDL_COUL && ir->fepvals->separate_dvdl[efptCOUL]) ||
243 (i == F_DVDL_VDW && ir->fepvals->separate_dvdl[efptVDW]) ||
244 (i == F_DVDL_BONDED && ir->fepvals->separate_dvdl[efptBONDED]) ||
245 (i == F_DVDL_RESTRAINT && ir->fepvals->separate_dvdl[efptRESTRAINT]) ||
246 (i == F_DKDL && ir->fepvals->separate_dvdl[efptMASS]) ||
247 (i == F_DVDL && ir->fepvals->separate_dvdl[efptFEP]))
249 md->bEner[i] = (ir->efep != efepNO);
251 else if ((interaction_function[i].flags & IF_VSITE) ||
252 (i == F_CONSTR) || (i == F_CONSTRNC) || (i == F_SETTLE))
254 md->bEner[i] = FALSE;
256 else if ((i == F_COUL_SR) || (i == F_EPOT) || (i == F_PRES) || (i == F_EQM))
260 else if ((i == F_GBPOL) && ir->implicit_solvent == eisGBSA)
264 else if ((i == F_NPSOLVATION) && ir->implicit_solvent == eisGBSA && (ir->sa_algorithm != esaNO))
268 else if ((i == F_GB12) || (i == F_GB13) || (i == F_GB14))
270 md->bEner[i] = FALSE;
272 else if ((i == F_ETOT) || (i == F_EKIN) || (i == F_TEMP))
274 md->bEner[i] = EI_DYNAMICS(ir->eI);
276 else if (i == F_DISPCORR || i == F_PDISPCORR)
278 md->bEner[i] = (ir->eDispCorr != edispcNO);
280 else if (i == F_DISRESVIOL)
282 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_DISRES) > 0);
284 else if (i == F_ORIRESDEV)
286 md->bEner[i] = (gmx_mtop_ftype_count(mtop, F_ORIRES) > 0);
288 else if (i == F_CONNBONDS)
290 md->bEner[i] = FALSE;
292 else if (i == F_COM_PULL)
294 md->bEner[i] = (ir->ePull == epullUMBRELLA || ir->ePull == epullCONST_F || ir->bRot);
296 else if (i == F_ECONSERVED)
298 md->bEner[i] = ((ir->etc == etcNOSEHOOVER || ir->etc == etcVRESCALE) &&
299 (ir->epc == epcNO || ir->epc == epcMTTK));
303 md->bEner[i] = (gmx_mtop_ftype_count(mtop, i) > 0);
307 /* OpenMM always produces only the following 4 energy terms */
308 md->bEner[F_EPOT] = TRUE;
309 md->bEner[F_EKIN] = TRUE;
310 md->bEner[F_ETOT] = TRUE;
311 md->bEner[F_TEMP] = TRUE;
314 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
315 if (ir->bAdress && !debug)
317 for (i = 0; i < F_NRE; i++)
319 md->bEner[i] = FALSE;
336 for (i = 0; i < F_NRE; i++)
340 ener_nm[md->f_nre] = interaction_function[i].longname;
346 md->bDiagPres = !TRICLINIC(ir->ref_p);
347 md->ref_p = (ir->ref_p[XX][XX]+ir->ref_p[YY][YY]+ir->ref_p[ZZ][ZZ])/DIM;
348 md->bTricl = TRICLINIC(ir->compress) || TRICLINIC(ir->deform);
349 md->bDynBox = DYNAMIC_BOX(*ir);
351 md->bNHC_trotter = IR_NVT_TROTTER(ir);
352 md->bPrintNHChains = ir->bPrintNHChains;
353 md->bMTTK = (IR_NPT_TROTTER(ir) || IR_NPH_TROTTER(ir));
354 md->bMu = NEED_MUTOT(*ir);
356 md->ebin = mk_ebin();
357 /* Pass NULL for unit to let get_ebin_space determine the units
358 * for interaction_function[i].longname
360 md->ie = get_ebin_space(md->ebin, md->f_nre, ener_nm, NULL);
363 /* This should be called directly after the call for md->ie,
364 * such that md->iconrmsd follows directly in the list.
366 md->iconrmsd = get_ebin_space(md->ebin, md->nCrmsd, conrmsd_nm, "");
370 md->ib = get_ebin_space(md->ebin,
371 md->bTricl ? NTRICLBOXS : NBOXS,
372 md->bTricl ? tricl_boxs_nm : boxs_nm,
374 md->ivol = get_ebin_space(md->ebin, 1, vol_nm, unit_volume);
375 md->idens = get_ebin_space(md->ebin, 1, dens_nm, unit_density_SI);
378 md->ipv = get_ebin_space(md->ebin, 1, pv_nm, unit_energy);
379 md->ienthalpy = get_ebin_space(md->ebin, 1, enthalpy_nm, unit_energy);
384 md->isvir = get_ebin_space(md->ebin, asize(sv_nm), sv_nm, unit_energy);
385 md->ifvir = get_ebin_space(md->ebin, asize(fv_nm), fv_nm, unit_energy);
389 md->ivir = get_ebin_space(md->ebin, asize(vir_nm), vir_nm, unit_energy);
393 md->ipres = get_ebin_space(md->ebin, asize(pres_nm), pres_nm, unit_pres_bar);
397 md->isurft = get_ebin_space(md->ebin, asize(surft_nm), surft_nm,
400 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
402 md->ipc = get_ebin_space(md->ebin, md->bTricl ? 6 : 3,
403 boxvel_nm, unit_vel);
407 md->imu = get_ebin_space(md->ebin, asize(mu_nm), mu_nm, unit_dipole_D);
409 if (ir->cos_accel != 0)
411 md->ivcos = get_ebin_space(md->ebin, asize(vcos_nm), vcos_nm, unit_vel);
412 md->ivisc = get_ebin_space(md->ebin, asize(visc_nm), visc_nm,
416 /* Energy monitoring */
417 for (i = 0; i < egNR; i++)
419 md->bEInd[i] = FALSE;
421 md->bEInd[egCOULSR] = TRUE;
422 md->bEInd[egLJSR ] = TRUE;
424 if (ir->rcoulomb > ir->rlist)
426 md->bEInd[egCOULLR] = TRUE;
430 if (ir->rvdw > ir->rlist)
432 md->bEInd[egLJLR] = TRUE;
437 md->bEInd[egLJSR] = FALSE;
438 md->bEInd[egBHAMSR] = TRUE;
439 if (ir->rvdw > ir->rlist)
441 md->bEInd[egBHAMLR] = TRUE;
446 md->bEInd[egLJ14] = TRUE;
447 md->bEInd[egCOUL14] = TRUE;
450 for (i = 0; (i < egNR); i++)
458 n = groups->grps[egcENER].nr;
459 /* for adress simulations, most energy terms are not meaningfull, and thus disabled*/
462 /*standard simulation*/
464 md->nE = (n*(n+1))/2;
468 /*AdResS simulation*/
475 snew(md->igrp, md->nE);
480 for (k = 0; (k < md->nEc); k++)
482 snew(gnm[k], STRLEN);
484 for (i = 0; (i < groups->grps[egcENER].nr); i++)
486 ni = groups->grps[egcENER].nm_ind[i];
487 for (j = i; (j < groups->grps[egcENER].nr); j++)
489 nj = groups->grps[egcENER].nm_ind[j];
490 for (k = kk = 0; (k < egNR); k++)
494 sprintf(gnm[kk], "%s:%s-%s", egrp_nm[k],
495 *(groups->grpname[ni]), *(groups->grpname[nj]));
499 md->igrp[n] = get_ebin_space(md->ebin, md->nEc,
500 (const char **)gnm, unit_energy);
504 for (k = 0; (k < md->nEc); k++)
512 gmx_incons("Number of energy terms wrong");
516 md->nTC = groups->grps[egcTC].nr;
517 md->nNHC = ir->opts.nhchainlength; /* shorthand for number of NH chains */
520 md->nTCP = 1; /* assume only one possible coupling system for barostat
527 if (md->etc == etcNOSEHOOVER)
529 if (md->bNHC_trotter)
531 md->mde_n = 2*md->nNHC*md->nTC;
535 md->mde_n = 2*md->nTC;
537 if (md->epc == epcMTTK)
539 md->mdeb_n = 2*md->nNHC*md->nTCP;
548 snew(md->tmp_r, md->mde_n);
549 snew(md->tmp_v, md->mde_n);
550 snew(md->grpnms, md->mde_n);
553 for (i = 0; (i < md->nTC); i++)
555 ni = groups->grps[egcTC].nm_ind[i];
556 sprintf(buf, "T-%s", *(groups->grpname[ni]));
557 grpnms[i] = strdup(buf);
559 md->itemp = get_ebin_space(md->ebin, md->nTC, (const char **)grpnms,
562 if (md->etc == etcNOSEHOOVER)
564 if (md->bPrintNHChains)
566 if (md->bNHC_trotter)
568 for (i = 0; (i < md->nTC); i++)
570 ni = groups->grps[egcTC].nm_ind[i];
571 bufi = *(groups->grpname[ni]);
572 for (j = 0; (j < md->nNHC); j++)
574 sprintf(buf, "Xi-%d-%s", j, bufi);
575 grpnms[2*(i*md->nNHC+j)] = strdup(buf);
576 sprintf(buf, "vXi-%d-%s", j, bufi);
577 grpnms[2*(i*md->nNHC+j)+1] = strdup(buf);
580 md->itc = get_ebin_space(md->ebin, md->mde_n,
581 (const char **)grpnms, unit_invtime);
584 for (i = 0; (i < md->nTCP); i++)
586 bufi = baro_nm[0]; /* All barostat DOF's together for now. */
587 for (j = 0; (j < md->nNHC); j++)
589 sprintf(buf, "Xi-%d-%s", j, bufi);
590 grpnms[2*(i*md->nNHC+j)] = strdup(buf);
591 sprintf(buf, "vXi-%d-%s", j, bufi);
592 grpnms[2*(i*md->nNHC+j)+1] = strdup(buf);
595 md->itcb = get_ebin_space(md->ebin, md->mdeb_n,
596 (const char **)grpnms, unit_invtime);
601 for (i = 0; (i < md->nTC); i++)
603 ni = groups->grps[egcTC].nm_ind[i];
604 bufi = *(groups->grpname[ni]);
605 sprintf(buf, "Xi-%s", bufi);
606 grpnms[2*i] = strdup(buf);
607 sprintf(buf, "vXi-%s", bufi);
608 grpnms[2*i+1] = strdup(buf);
610 md->itc = get_ebin_space(md->ebin, md->mde_n,
611 (const char **)grpnms, unit_invtime);
615 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
616 md->etc == etcVRESCALE)
618 for (i = 0; (i < md->nTC); i++)
620 ni = groups->grps[egcTC].nm_ind[i];
621 sprintf(buf, "Lamb-%s", *(groups->grpname[ni]));
622 grpnms[i] = strdup(buf);
624 md->itc = get_ebin_space(md->ebin, md->mde_n, (const char **)grpnms, "");
630 md->nU = groups->grps[egcACC].nr;
633 snew(grpnms, 3*md->nU);
634 for (i = 0; (i < md->nU); i++)
636 ni = groups->grps[egcACC].nm_ind[i];
637 sprintf(buf, "Ux-%s", *(groups->grpname[ni]));
638 grpnms[3*i+XX] = strdup(buf);
639 sprintf(buf, "Uy-%s", *(groups->grpname[ni]));
640 grpnms[3*i+YY] = strdup(buf);
641 sprintf(buf, "Uz-%s", *(groups->grpname[ni]));
642 grpnms[3*i+ZZ] = strdup(buf);
644 md->iu = get_ebin_space(md->ebin, 3*md->nU, (const char **)grpnms, unit_vel);
650 do_enxnms(fp_ene, &md->ebin->nener, &md->ebin->enm);
653 md->print_grpnms = NULL;
655 /* check whether we're going to write dh histograms */
657 if (ir->fepvals->separate_dhdl_file == esepdhdlfileNO)
659 /* Currently dh histograms are only written with dynamics */
660 if (EI_DYNAMICS(ir->eI))
664 mde_delta_h_coll_init(md->dhc, ir);
667 snew(md->dE, ir->fepvals->n_lambda);
671 md->fp_dhdl = fp_dhdl;
672 snew(md->dE, ir->fepvals->n_lambda);
677 snew(md->temperatures, ir->fepvals->n_lambda);
678 for (i = 0; i < ir->fepvals->n_lambda; i++)
680 md->temperatures[i] = ir->simtempvals->temperatures[i];
686 /* print a lambda vector to a string
687 fep = the inputrec's FEP input data
688 i = the index of the lambda vector
689 get_native_lambda = whether to print the native lambda
690 get_names = whether to print the names rather than the values
691 str = the pre-allocated string buffer to print to. */
692 static void print_lambda_vector(t_lambda *fep, int i,
693 gmx_bool get_native_lambda, gmx_bool get_names,
700 for (j = 0; j < efptNR; j++)
702 if (fep->separate_dvdl[j])
707 str[0] = 0; /* reset the string */
710 str += sprintf(str, "("); /* set the opening parenthesis*/
712 for (j = 0; j < efptNR; j++)
714 if (fep->separate_dvdl[j])
719 if (get_native_lambda && fep->init_lambda >= 0)
721 str += sprintf(str, "%.4f", fep->init_lambda);
725 str += sprintf(str, "%.4f", fep->all_lambda[j][i]);
730 str += sprintf(str, "%s", efpt_singular_names[j]);
732 /* print comma for the next item */
735 str += sprintf(str, ", ");
742 /* and add the closing parenthesis */
743 str += sprintf(str, ")");
748 extern FILE *open_dhdl(const char *filename, const t_inputrec *ir,
749 const output_env_t oenv)
752 const char *dhdl = "dH/d\\lambda", *deltag = "\\DeltaH", *lambda = "\\lambda",
753 *lambdastate = "\\lambda state", *remain = "remaining";
754 char title[STRLEN], label_x[STRLEN], label_y[STRLEN];
755 int i, np, nps, nsets, nsets_de, nsetsbegin;
756 int n_lambda_terms = 0;
757 t_lambda *fep = ir->fepvals; /* for simplicity */
758 t_expanded *expand = ir->expandedvals;
760 char buf[STRLEN], lambda_vec_str[STRLEN], lambda_name_str[STRLEN];
766 gmx_bool write_pV = FALSE;
768 /* count the number of different lambda terms */
769 for (i = 0; i < efptNR; i++)
771 if (fep->separate_dvdl[i])
777 if (fep->n_lambda == 0)
779 sprintf(title, "%s", dhdl);
780 sprintf(label_x, "Time (ps)");
781 sprintf(label_y, "%s (%s %s)",
782 dhdl, unit_energy, "[\\lambda]\\S-1\\N");
786 sprintf(title, "%s and %s", dhdl, deltag);
787 sprintf(label_x, "Time (ps)");
788 sprintf(label_y, "%s and %s (%s %s)",
789 dhdl, deltag, unit_energy, "[\\8l\\4]\\S-1\\N");
791 fp = gmx_fio_fopen(filename, "w+");
792 xvgr_header(fp, title, label_x, label_y, exvggtXNY, oenv);
796 bufplace = sprintf(buf, "T = %g (K) ",
799 if (ir->efep != efepSLOWGROWTH)
801 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
803 /* compatibility output */
804 sprintf(&(buf[bufplace]), "%s = %.4f", lambda, fep->init_lambda);
808 print_lambda_vector(fep, fep->init_fep_state, TRUE, FALSE,
810 print_lambda_vector(fep, fep->init_fep_state, TRUE, TRUE,
812 sprintf(&(buf[bufplace]), "%s %d: %s = %s",
813 lambdastate, fep->init_fep_state,
814 lambda_name_str, lambda_vec_str);
817 xvgr_subtitle(fp, buf, oenv);
821 if (fep->dhdl_derivatives == edhdlderivativesYES)
823 nsets_dhdl = n_lambda_terms;
825 /* count the number of delta_g states */
826 nsets_de = fep->lambda_stop_n - fep->lambda_start_n;
828 nsets = nsets_dhdl + nsets_de; /* dhdl + fep differences */
830 if (fep->n_lambda > 0 && (expand->elmcmove > elmcmoveNO))
832 nsets += 1; /*add fep state for expanded ensemble */
835 if (fep->bPrintEnergy)
837 nsets += 1; /* add energy to the dhdl as well */
841 if ((ir->epc != epcNO) && (fep->n_lambda > 0) && (fep->init_lambda < 0))
843 nsetsextend += 1; /* for PV term, other terms possible if required for
844 the reduced potential (only needed with foreign
845 lambda, and only output when init_lambda is not
846 set in order to maintain compatibility of the
850 snew(setname, nsetsextend);
852 if (expand->elmcmove > elmcmoveNO)
854 /* state for the fep_vals, if we have alchemical sampling */
855 sprintf(buf, "%s", "Thermodynamic state");
856 setname[s] = strdup(buf);
860 if (fep->bPrintEnergy)
862 sprintf(buf, "%s (%s)", "Energy", unit_energy);
863 setname[s] = strdup(buf);
867 if (fep->dhdl_derivatives == edhdlderivativesYES)
869 for (i = 0; i < efptNR; i++)
871 if (fep->separate_dvdl[i])
874 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
876 /* compatibility output */
877 sprintf(buf, "%s %s %.4f", dhdl, lambda, fep->init_lambda);
881 double lam = fep->init_lambda;
882 if (fep->init_lambda < 0)
884 lam = fep->all_lambda[i][fep->init_fep_state];
886 sprintf(buf, "%s %s = %.4f", dhdl, efpt_singular_names[i],
889 setname[s] = strdup(buf);
895 if (fep->n_lambda > 0)
897 /* g_bar has to determine the lambda values used in this simulation
898 * from this xvg legend.
901 if (expand->elmcmove > elmcmoveNO)
903 nsetsbegin = 1; /* for including the expanded ensemble */
910 if (fep->bPrintEnergy)
914 nsetsbegin += nsets_dhdl;
916 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
918 print_lambda_vector(fep, i, FALSE, FALSE, lambda_vec_str);
919 if ( (fep->init_lambda >= 0) && (n_lambda_terms == 1 ))
921 /* for compatible dhdl.xvg files */
922 nps = sprintf(buf, "%s %s %s", deltag, lambda, lambda_vec_str);
926 nps = sprintf(buf, "%s %s to %s", deltag, lambda, lambda_vec_str);
931 /* print the temperature for this state if doing simulated annealing */
932 sprintf(&buf[nps], "T = %g (%s)",
933 ir->simtempvals->temperatures[s-(nsetsbegin)],
936 setname[s] = strdup(buf);
941 np = sprintf(buf, "pV (%s)", unit_energy);
942 setname[nsetsextend-1] = strdup(buf); /* the first entry after
946 xvgr_legend(fp, nsetsextend, (const char **)setname, oenv);
948 for (s = 0; s < nsetsextend; s++)
958 static void copy_energy(t_mdebin *md, real e[], real ecpy[])
962 for (i = j = 0; (i < F_NRE); i++)
971 gmx_incons("Number of energy terms wrong");
975 void upd_mdebin(t_mdebin *md,
980 gmx_enerdata_t *enerd,
989 gmx_ekindata_t *ekind,
993 int i, j, k, kk, m, n, gid;
994 real crmsd[2], tmp6[6];
995 real bs[NTRICLBOXS], vol, dens, pv, enthalpy;
998 double store_dhdl[efptNR];
999 real store_energy = 0;
1002 /* Do NOT use the box in the state variable, but the separate box provided
1003 * as an argument. This is because we sometimes need to write the box from
1004 * the last timestep to match the trajectory frames.
1006 copy_energy(md, enerd->term, ecopy);
1007 add_ebin(md->ebin, md->ie, md->f_nre, ecopy, bSum);
1010 crmsd[0] = constr_rmsd(constr, FALSE);
1013 crmsd[1] = constr_rmsd(constr, TRUE);
1015 add_ebin(md->ebin, md->iconrmsd, md->nCrmsd, crmsd, FALSE);
1022 bs[0] = box[XX][XX];
1023 bs[1] = box[YY][YY];
1024 bs[2] = box[ZZ][ZZ];
1025 bs[3] = box[YY][XX];
1026 bs[4] = box[ZZ][XX];
1027 bs[5] = box[ZZ][YY];
1032 bs[0] = box[XX][XX];
1033 bs[1] = box[YY][YY];
1034 bs[2] = box[ZZ][ZZ];
1037 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1038 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1039 add_ebin(md->ebin, md->ib, nboxs, bs, bSum);
1040 add_ebin(md->ebin, md->ivol, 1, &vol, bSum);
1041 add_ebin(md->ebin, md->idens, 1, &dens, bSum);
1045 /* This is pV (in kJ/mol). The pressure is the reference pressure,
1046 not the instantaneous pressure */
1047 pv = vol*md->ref_p/PRESFAC;
1049 add_ebin(md->ebin, md->ipv, 1, &pv, bSum);
1050 enthalpy = pv + enerd->term[F_ETOT];
1051 add_ebin(md->ebin, md->ienthalpy, 1, &enthalpy, bSum);
1056 add_ebin(md->ebin, md->isvir, 9, svir[0], bSum);
1057 add_ebin(md->ebin, md->ifvir, 9, fvir[0], bSum);
1061 add_ebin(md->ebin, md->ivir, 9, vir[0], bSum);
1065 add_ebin(md->ebin, md->ipres, 9, pres[0], bSum);
1069 tmp = (pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ];
1070 add_ebin(md->ebin, md->isurft, 1, &tmp, bSum);
1072 if (md->epc == epcPARRINELLORAHMAN || md->epc == epcMTTK)
1074 tmp6[0] = state->boxv[XX][XX];
1075 tmp6[1] = state->boxv[YY][YY];
1076 tmp6[2] = state->boxv[ZZ][ZZ];
1077 tmp6[3] = state->boxv[YY][XX];
1078 tmp6[4] = state->boxv[ZZ][XX];
1079 tmp6[5] = state->boxv[ZZ][YY];
1080 add_ebin(md->ebin, md->ipc, md->bTricl ? 6 : 3, tmp6, bSum);
1084 add_ebin(md->ebin, md->imu, 3, mu_tot, bSum);
1086 if (ekind && ekind->cosacc.cos_accel != 0)
1088 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
1089 dens = (tmass*AMU)/(vol*NANO*NANO*NANO);
1090 add_ebin(md->ebin, md->ivcos, 1, &(ekind->cosacc.vcos), bSum);
1091 /* 1/viscosity, unit 1/(kg m^-1 s^-1) */
1092 tmp = 1/(ekind->cosacc.cos_accel/(ekind->cosacc.vcos*PICO)
1093 *dens*vol*sqr(box[ZZ][ZZ]*NANO/(2*M_PI)));
1094 add_ebin(md->ebin, md->ivisc, 1, &tmp, bSum);
1099 for (i = 0; (i < md->nEg); i++)
1101 for (j = i; (j < md->nEg); j++)
1103 gid = GID(i, j, md->nEg);
1104 for (k = kk = 0; (k < egNR); k++)
1108 eee[kk++] = enerd->grpp.ener[k][gid];
1111 add_ebin(md->ebin, md->igrp[n], md->nEc, eee, bSum);
1119 for (i = 0; (i < md->nTC); i++)
1121 md->tmp_r[i] = ekind->tcstat[i].T;
1123 add_ebin(md->ebin, md->itemp, md->nTC, md->tmp_r, bSum);
1125 if (md->etc == etcNOSEHOOVER)
1127 /* whether to print Nose-Hoover chains: */
1128 if (md->bPrintNHChains)
1130 if (md->bNHC_trotter)
1132 for (i = 0; (i < md->nTC); i++)
1134 for (j = 0; j < md->nNHC; j++)
1137 md->tmp_r[2*k] = state->nosehoover_xi[k];
1138 md->tmp_r[2*k+1] = state->nosehoover_vxi[k];
1141 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1145 for (i = 0; (i < md->nTCP); i++)
1147 for (j = 0; j < md->nNHC; j++)
1150 md->tmp_r[2*k] = state->nhpres_xi[k];
1151 md->tmp_r[2*k+1] = state->nhpres_vxi[k];
1154 add_ebin(md->ebin, md->itcb, md->mdeb_n, md->tmp_r, bSum);
1159 for (i = 0; (i < md->nTC); i++)
1161 md->tmp_r[2*i] = state->nosehoover_xi[i];
1162 md->tmp_r[2*i+1] = state->nosehoover_vxi[i];
1164 add_ebin(md->ebin, md->itc, md->mde_n, md->tmp_r, bSum);
1168 else if (md->etc == etcBERENDSEN || md->etc == etcYES ||
1169 md->etc == etcVRESCALE)
1171 for (i = 0; (i < md->nTC); i++)
1173 md->tmp_r[i] = ekind->tcstat[i].lambda;
1175 add_ebin(md->ebin, md->itc, md->nTC, md->tmp_r, bSum);
1179 if (ekind && md->nU > 1)
1181 for (i = 0; (i < md->nU); i++)
1183 copy_rvec(ekind->grpstat[i].u, md->tmp_v[i]);
1185 add_ebin(md->ebin, md->iu, 3*md->nU, md->tmp_v[0], bSum);
1188 ebin_increase_count(md->ebin, bSum);
1190 /* BAR + thermodynamic integration values */
1191 if ((md->fp_dhdl || md->dhc) && bDoDHDL)
1193 for (i = 0; i < enerd->n_lambda-1; i++)
1195 /* zero for simulated tempering */
1196 md->dE[i] = enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0];
1197 if (md->temperatures != NULL)
1199 /* MRS: is this right, given the way we have defined the exchange probabilities? */
1200 /* is this even useful to have at all? */
1201 md->dE[i] += (md->temperatures[i]/
1202 md->temperatures[state->fep_state]-1.0)*
1203 enerd->term[F_EKIN];
1209 fprintf(md->fp_dhdl, "%.4f", time);
1210 /* the current free energy state */
1212 /* print the current state if we are doing expanded ensemble */
1213 if (expand->elmcmove > elmcmoveNO)
1215 fprintf(md->fp_dhdl, " %4d", state->fep_state);
1217 /* total energy (for if the temperature changes */
1218 if (fep->bPrintEnergy)
1220 store_energy = enerd->term[F_ETOT];
1221 fprintf(md->fp_dhdl, " %#.8g", store_energy);
1224 if (fep->dhdl_derivatives == edhdlderivativesYES)
1226 for (i = 0; i < efptNR; i++)
1228 if (fep->separate_dvdl[i])
1230 /* assumes F_DVDL is first */
1231 fprintf(md->fp_dhdl, " %#.8g", enerd->term[F_DVDL+i]);
1235 for (i = fep->lambda_start_n; i < fep->lambda_stop_n; i++)
1237 fprintf(md->fp_dhdl, " %#.8g", md->dE[i]);
1239 if ((md->epc != epcNO) &&
1240 (enerd->n_lambda > 0) &&
1241 (fep->init_lambda < 0))
1243 fprintf(md->fp_dhdl, " %#.8g", pv); /* PV term only needed when
1244 there are alternate state
1245 lambda and we're not in
1246 compatibility mode */
1248 fprintf(md->fp_dhdl, "\n");
1249 /* and the binary free energy output */
1251 if (md->dhc && bDoDHDL)
1254 for (i = 0; i < efptNR; i++)
1256 if (fep->separate_dvdl[i])
1258 /* assumes F_DVDL is first */
1259 store_dhdl[idhdl] = enerd->term[F_DVDL+i];
1263 store_energy = enerd->term[F_ETOT];
1264 /* store_dh is dE */
1265 mde_delta_h_coll_add_dh(md->dhc,
1266 (double)state->fep_state,
1270 md->dE + fep->lambda_start_n,
1277 void upd_mdebin_step(t_mdebin *md)
1279 ebin_increase_count(md->ebin, FALSE);
1282 static void npr(FILE *log, int n, char c)
1284 for (; (n > 0); n--)
1286 fprintf(log, "%c", c);
1290 static void pprint(FILE *log, const char *s, t_mdebin *md)
1294 char buf1[22], buf2[22];
1297 fprintf(log, "\t<====== ");
1298 npr(log, slen, CHAR);
1299 fprintf(log, " ==>\n");
1300 fprintf(log, "\t<==== %s ====>\n", s);
1301 fprintf(log, "\t<== ");
1302 npr(log, slen, CHAR);
1303 fprintf(log, " ======>\n\n");
1305 fprintf(log, "\tStatistics over %s steps using %s frames\n",
1306 gmx_step_str(md->ebin->nsteps_sim, buf1),
1307 gmx_step_str(md->ebin->nsum_sim, buf2));
1311 void print_ebin_header(FILE *log, gmx_large_int_t steps, double time, real lambda)
1315 fprintf(log, " %12s %12s %12s\n"
1316 " %12s %12.5f %12.5f\n\n",
1317 "Step", "Time", "Lambda", gmx_step_str(steps, buf), time, lambda);
1320 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
1322 gmx_large_int_t step, double time,
1323 int mode, gmx_bool bCompact,
1324 t_mdebin *md, t_fcdata *fcd,
1325 gmx_groups_t *groups, t_grpopts *opts)
1327 /*static char **grpnms=NULL;*/
1329 int i, j, n, ni, nj, ndr, nor, b;
1331 real *disre_rm3tav, *disre_rt;
1333 /* these are for the old-style blocks (1 subblock, only reals), because
1334 there can be only one per ID for these */
1339 /* temporary arrays for the lambda values to write out */
1340 double enxlambda_data[2];
1350 fr.nsteps = md->ebin->nsteps;
1351 fr.dt = md->delta_t;
1352 fr.nsum = md->ebin->nsum;
1353 fr.nre = (bEne) ? md->ebin->nener : 0;
1354 fr.ener = md->ebin->e;
1355 ndisre = bDR ? fcd->disres.npair : 0;
1356 disre_rm3tav = fcd->disres.rm3tav;
1357 disre_rt = fcd->disres.rt;
1358 /* Optional additional old-style (real-only) blocks. */
1359 for (i = 0; i < enxNR; i++)
1363 if (fcd->orires.nr > 0 && bOR)
1365 diagonalize_orires_tensors(&(fcd->orires));
1366 nr[enxOR] = fcd->orires.nr;
1367 block[enxOR] = fcd->orires.otav;
1369 nr[enxORI] = (fcd->orires.oinsl != fcd->orires.otav) ?
1371 block[enxORI] = fcd->orires.oinsl;
1372 id[enxORI] = enxORI;
1373 nr[enxORT] = fcd->orires.nex*12;
1374 block[enxORT] = fcd->orires.eig;
1375 id[enxORT] = enxORT;
1378 /* whether we are going to wrte anything out: */
1379 if (fr.nre || ndisre || nr[enxOR] || nr[enxORI])
1382 /* the old-style blocks go first */
1384 for (i = 0; i < enxNR; i++)
1391 add_blocks_enxframe(&fr, fr.nblock);
1392 for (b = 0; b < fr.nblock; b++)
1394 add_subblocks_enxblock(&(fr.block[b]), 1);
1395 fr.block[b].id = id[b];
1396 fr.block[b].sub[0].nr = nr[b];
1398 fr.block[b].sub[0].type = xdr_datatype_float;
1399 fr.block[b].sub[0].fval = block[b];
1401 fr.block[b].sub[0].type = xdr_datatype_double;
1402 fr.block[b].sub[0].dval = block[b];
1406 /* check for disre block & fill it. */
1411 add_blocks_enxframe(&fr, fr.nblock);
1413 add_subblocks_enxblock(&(fr.block[db]), 2);
1414 fr.block[db].id = enxDISRE;
1415 fr.block[db].sub[0].nr = ndisre;
1416 fr.block[db].sub[1].nr = ndisre;
1418 fr.block[db].sub[0].type = xdr_datatype_float;
1419 fr.block[db].sub[1].type = xdr_datatype_float;
1420 fr.block[db].sub[0].fval = disre_rt;
1421 fr.block[db].sub[1].fval = disre_rm3tav;
1423 fr.block[db].sub[0].type = xdr_datatype_double;
1424 fr.block[db].sub[1].type = xdr_datatype_double;
1425 fr.block[db].sub[0].dval = disre_rt;
1426 fr.block[db].sub[1].dval = disre_rm3tav;
1429 /* here we can put new-style blocks */
1431 /* Free energy perturbation blocks */
1434 mde_delta_h_coll_handle_block(md->dhc, &fr, fr.nblock);
1437 /* we can now free & reset the data in the blocks */
1440 mde_delta_h_coll_reset(md->dhc);
1443 /* do the actual I/O */
1444 do_enx(fp_ene, &fr);
1445 gmx_fio_check_file_position(enx_file_pointer(fp_ene));
1448 /* We have stored the sums, so reset the sum history */
1449 reset_ebin_sums(md->ebin);
1457 pprint(log, "A V E R A G E S", md);
1463 pprint(log, "R M S - F L U C T U A T I O N S", md);
1467 gmx_fatal(FARGS, "Invalid print mode (%d)", mode);
1472 for (i = 0; i < opts->ngtc; i++)
1474 if (opts->annealing[i] != eannNO)
1476 fprintf(log, "Current ref_t for group %s: %8.1f\n",
1477 *(groups->grpname[groups->grps[egcTC].nm_ind[i]]),
1481 if (mode == eprNORMAL && fcd->orires.nr > 0)
1483 print_orires_log(log, &(fcd->orires));
1485 fprintf(log, " Energies (%s)\n", unit_energy);
1486 pr_ebin(log, md->ebin, md->ie, md->f_nre+md->nCrmsd, 5, mode, TRUE);
1493 pr_ebin(log, md->ebin, md->ib, md->bTricl ? NTRICLBOXS : NBOXS, 5,
1499 fprintf(log, " Constraint Virial (%s)\n", unit_energy);
1500 pr_ebin(log, md->ebin, md->isvir, 9, 3, mode, FALSE);
1502 fprintf(log, " Force Virial (%s)\n", unit_energy);
1503 pr_ebin(log, md->ebin, md->ifvir, 9, 3, mode, FALSE);
1508 fprintf(log, " Total Virial (%s)\n", unit_energy);
1509 pr_ebin(log, md->ebin, md->ivir, 9, 3, mode, FALSE);
1514 fprintf(log, " Pressure (%s)\n", unit_pres_bar);
1515 pr_ebin(log, md->ebin, md->ipres, 9, 3, mode, FALSE);
1520 fprintf(log, " Total Dipole (%s)\n", unit_dipole_D);
1521 pr_ebin(log, md->ebin, md->imu, 3, 3, mode, FALSE);
1527 if (md->print_grpnms == NULL)
1529 snew(md->print_grpnms, md->nE);
1531 for (i = 0; (i < md->nEg); i++)
1533 ni = groups->grps[egcENER].nm_ind[i];
1534 for (j = i; (j < md->nEg); j++)
1536 nj = groups->grps[egcENER].nm_ind[j];
1537 sprintf(buf, "%s-%s", *(groups->grpname[ni]),
1538 *(groups->grpname[nj]));
1539 md->print_grpnms[n++] = strdup(buf);
1543 sprintf(buf, "Epot (%s)", unit_energy);
1544 fprintf(log, "%15s ", buf);
1545 for (i = 0; (i < egNR); i++)
1549 fprintf(log, "%12s ", egrp_nm[i]);
1553 for (i = 0; (i < md->nE); i++)
1555 fprintf(log, "%15s", md->print_grpnms[i]);
1556 pr_ebin(log, md->ebin, md->igrp[i], md->nEc, md->nEc, mode,
1563 pr_ebin(log, md->ebin, md->itemp, md->nTC, 4, mode, TRUE);
1568 fprintf(log, "%15s %12s %12s %12s\n",
1569 "Group", "Ux", "Uy", "Uz");
1570 for (i = 0; (i < md->nU); i++)
1572 ni = groups->grps[egcACC].nm_ind[i];
1573 fprintf(log, "%15s", *groups->grpname[ni]);
1574 pr_ebin(log, md->ebin, md->iu+3*i, 3, 3, mode, FALSE);
1583 void update_energyhistory(energyhistory_t * enerhist, t_mdebin * mdebin)
1587 enerhist->nsteps = mdebin->ebin->nsteps;
1588 enerhist->nsum = mdebin->ebin->nsum;
1589 enerhist->nsteps_sim = mdebin->ebin->nsteps_sim;
1590 enerhist->nsum_sim = mdebin->ebin->nsum_sim;
1591 enerhist->nener = mdebin->ebin->nener;
1593 if (mdebin->ebin->nsum > 0)
1595 /* Check if we need to allocate first */
1596 if (enerhist->ener_ave == NULL)
1598 snew(enerhist->ener_ave, enerhist->nener);
1599 snew(enerhist->ener_sum, enerhist->nener);
1602 for (i = 0; i < enerhist->nener; i++)
1604 enerhist->ener_ave[i] = mdebin->ebin->e[i].eav;
1605 enerhist->ener_sum[i] = mdebin->ebin->e[i].esum;
1609 if (mdebin->ebin->nsum_sim > 0)
1611 /* Check if we need to allocate first */
1612 if (enerhist->ener_sum_sim == NULL)
1614 snew(enerhist->ener_sum_sim, enerhist->nener);
1617 for (i = 0; i < enerhist->nener; i++)
1619 enerhist->ener_sum_sim[i] = mdebin->ebin->e_sim[i].esum;
1624 mde_delta_h_coll_update_energyhistory(mdebin->dhc, enerhist);
1628 void restore_energyhistory_from_state(t_mdebin * mdebin,
1629 energyhistory_t * enerhist)
1633 if ((enerhist->nsum > 0 || enerhist->nsum_sim > 0) &&
1634 mdebin->ebin->nener != enerhist->nener)
1636 gmx_fatal(FARGS, "Mismatch between number of energies in run input (%d) and checkpoint file (%d).",
1637 mdebin->ebin->nener, enerhist->nener);
1640 mdebin->ebin->nsteps = enerhist->nsteps;
1641 mdebin->ebin->nsum = enerhist->nsum;
1642 mdebin->ebin->nsteps_sim = enerhist->nsteps_sim;
1643 mdebin->ebin->nsum_sim = enerhist->nsum_sim;
1645 for (i = 0; i < mdebin->ebin->nener; i++)
1647 mdebin->ebin->e[i].eav =
1648 (enerhist->nsum > 0 ? enerhist->ener_ave[i] : 0);
1649 mdebin->ebin->e[i].esum =
1650 (enerhist->nsum > 0 ? enerhist->ener_sum[i] : 0);
1651 mdebin->ebin->e_sim[i].esum =
1652 (enerhist->nsum_sim > 0 ? enerhist->ener_sum_sim[i] : 0);
1656 mde_delta_h_coll_restore_energyhistory(mdebin->dhc, enerhist);