-/* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
+/*
+ * This file is part of the GROMACS molecular simulation package.
*
+ * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
*
- * This source code is part of
- *
- * G R O M A C S
- *
- * GROningen MAchine for Chemical Simulations
- *
- * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
- * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
- * Copyright (c) 2001-2012, The GROMACS development team,
- * check out http://www.gromacs.org for more information.
- *
- * This program is free software; you can redistribute it and/or
- * modify it under the terms of the GNU General Public License
- * as published by the Free Software Foundation; either version 2
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
* of the License, or (at your option) any later version.
*
- * If you want to redistribute modifications, please consider that
- * scientific software is very special. Version control is crucial -
- * bugs must be traceable. We will be happy to consider code for
- * inclusion in the official distribution, but derived work must not
- * be called official GROMACS. Details are found in the README & COPYING
- * files - if they are missing, get the official version at www.gromacs.org.
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
+ * Lesser General Public License for more details.
*
- * To help us fund GROMACS development, we humbly ask that you cite
- * the papers on the package - you can find them in the top README file.
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
*
- * For more info, check our website at http://www.gromacs.org
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
*
- * And Hey:
- * GROwing Monsters And Cloning Shrimps
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
*/
-#ifdef HAVE_CONFIG_H
-#include <config.h>
-#endif
-
-#ifdef GMX_CRAY_XT3
-#include <catamount/dclock.h>
-#endif
+#include "gmxpre.h"
-
-#include <stdio.h>
-#include <time.h>
-#ifdef HAVE_SYS_TIME_H
-#include <sys/time.h>
-#endif
#include <math.h>
-#include "typedefs.h"
-#include "string2.h"
-#include "gmxfio.h"
-#include "smalloc.h"
-#include "names.h"
-#include "confio.h"
-#include "mvdata.h"
-#include "txtdump.h"
-#include "pbc.h"
-#include "chargegroup.h"
-#include "vec.h"
-#include "nrnb.h"
-#include "mshift.h"
-#include "mdrun.h"
-#include "update.h"
-#include "physics.h"
-#include "main.h"
-#include "mdatoms.h"
-#include "force.h"
-#include "bondf.h"
-#include "pme.h"
-#include "disre.h"
-#include "orires.h"
-#include "network.h"
-#include "calcmu.h"
-#include "constr.h"
-#include "xvgr.h"
-#include "trnio.h"
-#include "xtcio.h"
-#include "gmx_random.h"
-#include "domdec.h"
-#include "partdec.h"
-#include "gmx_wallcycle.h"
-#include "macros.h"
+#include <stdio.h>
+#include "gromacs/fileio/confio.h"
+#include "gromacs/fileio/gmxfio.h"
+#include "gromacs/fileio/trnio.h"
+#include "gromacs/fileio/xtcio.h"
+#include "gromacs/legacyheaders/calcmu.h"
+#include "gromacs/legacyheaders/chargegroup.h"
+#include "gromacs/legacyheaders/constr.h"
+#include "gromacs/legacyheaders/disre.h"
+#include "gromacs/legacyheaders/domdec.h"
+#include "gromacs/legacyheaders/force.h"
+#include "gromacs/legacyheaders/macros.h"
+#include "gromacs/legacyheaders/mdatoms.h"
+#include "gromacs/legacyheaders/mdrun.h"
+#include "gromacs/legacyheaders/names.h"
+#include "gromacs/legacyheaders/network.h"
+#include "gromacs/legacyheaders/nrnb.h"
+#include "gromacs/legacyheaders/orires.h"
+#include "gromacs/legacyheaders/pme.h"
+#include "gromacs/legacyheaders/txtdump.h"
+#include "gromacs/legacyheaders/typedefs.h"
+#include "gromacs/legacyheaders/update.h"
+#include "gromacs/math/units.h"
+#include "gromacs/math/vec.h"
+#include "gromacs/random/random.h"
+#include "gromacs/timing/wallcycle.h"
+#include "gromacs/utility/fatalerror.h"
#include "gromacs/utility/gmxmpi.h"
+#include "gromacs/utility/smalloc.h"
+
+static void init_df_history_weights(df_history_t *dfhist, t_expanded *expand, int nlim)
+{
+ int i;
+ dfhist->wl_delta = expand->init_wl_delta;
+ for (i = 0; i < nlim; i++)
+ {
+ dfhist->sum_weights[i] = expand->init_lambda_weights[i];
+ dfhist->sum_dg[i] = expand->init_lambda_weights[i];
+ }
+}
+
+/* Eventually should contain all the functions needed to initialize expanded ensemble
+ before the md loop starts */
+extern void init_expanded_ensemble(gmx_bool bStateFromCP, t_inputrec *ir, df_history_t *dfhist)
+{
+ if (!bStateFromCP)
+ {
+ init_df_history_weights(dfhist, ir->expandedvals, ir->fepvals->n_lambda);
+ }
+}
-void GenerateGibbsProbabilities(real *ene, real *p_k, real *pks, int minfep, int maxfep)
+static void GenerateGibbsProbabilities(real *ene, double *p_k, double *pks, int minfep, int maxfep)
{
int i;
}
}
-void GenerateWeightedGibbsProbabilities(real *ene, real *p_k, real *pks, int nlim, real *nvals, real delta)
+static void GenerateWeightedGibbsProbabilities(real *ene, double *p_k, double *pks, int nlim, real *nvals, real delta)
{
int i;
return bIfFlat;
}
-static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_large_int_t step)
+static gmx_bool CheckIfDoneEquilibrating(int nlim, t_expanded *expand, df_history_t *dfhist, gmx_int64_t step)
{
int i, totalsamples;
}
static gmx_bool UpdateWeights(int nlim, t_expanded *expand, df_history_t *dfhist,
- int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_large_int_t step)
+ int fep_state, real *scaled_lamee, real *weighted_lamee, gmx_int64_t step)
{
real maxdiff = 0.000000001;
gmx_bool bSufficientSamples;
- int i, k, n, nz, indexi, indexk, min_n, max_n, nlam, totali;
+ int i, k, n, nz, indexi, indexk, min_n, max_n, totali;
int n0, np1, nm1, nval, min_nvalm, min_nvalp, maxc;
real omega_m1_0, omega_p1_m1, omega_m1_p1, omega_p1_0, clam_osum;
- real de, de_function, dr, denom, maxdr, pks = 0;
+ real de, de_function, dr, denom, maxdr;
real min_val, cnval, zero_sum_weights;
real *omegam_array, *weightsm_array, *omegap_array, *weightsp_array, *varm_array, *varp_array, *dwp_array, *dwm_array;
real clam_varm, clam_varp, clam_weightsm, clam_weightsp, clam_minvar;
- real *lam_weights, *lam_minvar_corr, *lam_variance, *lam_dg, *p_k;
+ real *lam_weights, *lam_minvar_corr, *lam_variance, *lam_dg;
+ double *p_k;
+ double pks = 0;
real *numweighted_lamee, *logfrac;
int *nonzero;
real chi_m1_0, chi_p1_0, chi_m2_0, chi_p2_0, chi_p1_m1, chi_p2_m1, chi_m1_p1, chi_m2_p1;
GenerateGibbsProbabilities(weighted_lamee, p_k, &pks, 0, nlim-1);
for (i = 0; i < nlim; i++)
{
- dfhist->wl_histo[i] += p_k[i];
+ dfhist->wl_histo[i] += (real)p_k[i];
}
/* then increment weights (uses count) */
for (i = 0; i < nlim; i++)
{
- dfhist->sum_weights[i] -= dfhist->wl_delta*p_k[i];
+ dfhist->sum_weights[i] -= dfhist->wl_delta*(real)p_k[i];
}
/* Alternate definition, using logarithms. Shouldn't make very much difference! */
/*
real di;
for (i=0;i<nlim;i++)
{
- di = 1+dfhist->wl_delta*p_k[i];
+ di = (real)1.0 + dfhist->wl_delta*(real)p_k[i];
dfhist->sum_weights[i] -= log(di);
}
*/
return FALSE;
}
-static int ChooseNewLambda(int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, real *p_k, gmx_rng_t rng)
+static int ChooseNewLambda(int nlim, t_expanded *expand, df_history_t *dfhist, int fep_state, real *weighted_lamee, double *p_k,
+ gmx_int64_t seed, gmx_int64_t step)
{
/* Choose new lambda value, and update transition matrix */
int i, ifep, jfep, minfep, maxfep, lamnew, lamtrial, starting_fep_state;
- real r1, r2, pks, de_old, de_new, de, trialprob, tprob = 0;
+ real r1, r2, de_old, de_new, de, trialprob, tprob = 0;
real **Tij;
- real *propose, *accept, *remainder;
+ double *propose, *accept, *remainder;
+ double pks;
real sum, pnorm;
gmx_bool bRestricted;
for (i = 0; i < expand->lmc_repeats; i++)
{
+ double rnd[2];
+
+ gmx_rng_cycle_2uniform(step, i, seed, RND_SEED_EXPANDED, rnd);
for (ifep = 0; ifep < nlim; ifep++)
{
accept[ifep] = 1.0;
}
/* Gibbs sampling */
- r1 = gmx_rng_uniform_real(rng);
+ r1 = rnd[0];
for (lamnew = minfep; lamnew <= maxfep; lamnew++)
{
if (r1 <= p_k[lamnew])
}
}
- r1 = gmx_rng_uniform_real(rng);
+ r1 = rnd[0];
for (lamtrial = minfep; lamtrial <= maxfep; lamtrial++)
{
pnorm = p_k[lamtrial]/remainder[fep_state];
{
tprob = trialprob;
}
- r2 = gmx_rng_uniform_real(rng);
+ r2 = rnd[1];
if (r2 < tprob)
{
lamnew = lamtrial;
if (lamnew > maxfep)
{
/* it's possible some rounding is failing */
- if (remainder[fep_state] < 2.0e-15)
+ if (gmx_within_tol(remainder[fep_state], 0, 50*GMX_DOUBLE_EPS))
{
- /* probably numerical rounding error -- no state other than the original has weight */
+ /* numerical rounding error -- no state other than the original has weight */
lamnew = fep_state;
}
else
else if ((expand->elmcmove == elmcmoveMETROPOLIS) || (expand->elmcmove == elmcmoveBARKER))
{
/* use the metropolis sampler with trial +/- 1 */
- r1 = gmx_rng_uniform_real(rng);
+ r1 = rnd[0];
if (r1 < 0.5)
{
if (fep_state == 0)
accept[lamtrial] = 1.0;
}
- r2 = gmx_rng_uniform_real(rng);
+ r2 = rnd[1];
if (r2 < tprob)
{
lamnew = lamtrial;
/* print out the weights to the log, along with current state */
extern void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
- int nlam, int frequency, gmx_large_int_t step)
+ int fep_state, int frequency, gmx_int64_t step)
{
int nlim, i, ifep, jfep;
real dw, dg, dv, dm, Tprint;
{
fprintf(outfile, " %10.5f %10.5f", dfhist->sum_weights[ifep], dw);
}
- if (ifep == nlam)
+ if (ifep == fep_state)
{
fprintf(outfile, " <<\n");
}
}
}
-extern void get_mc_state(gmx_rng_t rng, t_state *state)
-{
- gmx_rng_get_state(rng, state->mc_rng, state->mc_rngi);
-}
-
-extern void set_mc_state(gmx_rng_t rng, t_state *state)
-{
- gmx_rng_set_state(rng, state->mc_rng, state->mc_rngi[0]);
-}
-
extern int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
- t_state *state, t_extmass *MassQ, df_history_t *dfhist,
- gmx_large_int_t step, gmx_rng_t mcrng,
+ t_state *state, t_extmass *MassQ, int fep_state, df_history_t *dfhist,
+ gmx_int64_t step,
rvec *v, t_mdatoms *mdatoms)
+/* Note that the state variable is only needed for simulated tempering, not
+ Hamiltonian expanded ensemble. May be able to remove it after integrator refactoring. */
{
- real *pfep_lamee, *p_k, *scaled_lamee, *weighted_lamee;
- int i, nlam, nlim, lamnew, totalsamples;
+ real *pfep_lamee, *scaled_lamee, *weighted_lamee;
+ double *p_k;
+ int i, nlim, lamnew, totalsamples;
real oneovert, maxscaled = 0, maxweighted = 0;
t_expanded *expand;
t_simtemp *simtemp;
expand = ir->expandedvals;
simtemp = ir->simtempvals;
nlim = ir->fepvals->n_lambda;
- nlam = state->fep_state;
snew(scaled_lamee, nlim);
snew(weighted_lamee, nlim);
snew(pfep_lamee, nlim);
snew(p_k, nlim);
- if (expand->bInit_weights) /* if initialized weights, we need to fill them in */
- {
- dfhist->wl_delta = expand->init_wl_delta; /* MRS -- this would fit better somewhere else? */
- for (i = 0; i < nlim; i++)
- {
- dfhist->sum_weights[i] = expand->init_lambda_weights[i];
- dfhist->sum_dg[i] = expand->init_lambda_weights[i];
- }
- expand->bInit_weights = FALSE;
- }
-
/* update the count at the current lambda*/
- dfhist->n_at_lam[nlam]++;
+ dfhist->n_at_lam[fep_state]++;
/* need to calculate the PV term somewhere, but not needed here? Not until there's a lambda state that's
pressure controlled.*/
{
/* Note -- this assumes no mass changes, since kinetic energy is not added . . . */
scaled_lamee[i] = (enerd->enerpart_lambda[i+1]-enerd->enerpart_lambda[0])/(simtemp->temperatures[i]*BOLTZ)
- + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[nlam]))/BOLTZ;
+ + enerd->term[F_EPOT]*(1.0/(simtemp->temperatures[i])- 1.0/(simtemp->temperatures[fep_state]))/BOLTZ;
}
else
{
{
for (i = 0; i < nlim; i++)
{
- scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[nlam])/BOLTZ;
+ scaled_lamee[i] = enerd->term[F_EPOT]*(1.0/simtemp->temperatures[i] - 1.0/simtemp->temperatures[fep_state])/BOLTZ;
}
}
}
/* update weights - we decide whether or not to actually do this inside */
- bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, nlam, scaled_lamee, weighted_lamee, step);
+ bDoneEquilibrating = UpdateWeights(nlim, expand, dfhist, fep_state, scaled_lamee, weighted_lamee, step);
if (bDoneEquilibrating)
{
if (log)
}
}
- lamnew = ChooseNewLambda(nlim, expand, dfhist, nlam, weighted_lamee, p_k, mcrng);
+ lamnew = ChooseNewLambda(nlim, expand, dfhist, fep_state, weighted_lamee, p_k,
+ ir->expandedvals->lmc_seed, step);
/* if using simulated tempering, we need to adjust the temperatures */
- if (ir->bSimTemp && (lamnew != nlam)) /* only need to change the temperatures if we change the state */
+ if (ir->bSimTemp && (lamnew != fep_state)) /* only need to change the temperatures if we change the state */
{
int i, j, n, d;
real *buf_ngtc;
/* we don't need to manipulate the ekind information, as it isn't due to be reset until the next step anyway */
- nstart = mdatoms->start;
- nend = nstart + mdatoms->homenr;
+ nstart = 0;
+ nend = mdatoms->homenr;
for (n = nstart; n < nend; n++)
{
gt = 0;