WHAM package
Submodules
WHAM.binless module
- Implementation of binless WHAM (described in Shirts M., & Chodera J. D. (2008)) using
Negative log-likelihood maximization, inspired by Zhu, F., & Hummer, G. (2012).
Self-consistent iteration.
- class WHAM.binless.Calc1D
Bases:
WHAM.binless.CalcBaseClass containing methods to compute free energy profiles from 1D umbrella sampling data (i.e. data from biasing a single order parameter) using binless WHAM.
- x_l
1-dimensional array containing unrolled order parameter which is being biased.
- Type
ndarray
- G_l
1-dimensional array containing WHAM-computed weights corresponding to each (unrolled) order parameter.
- Type
ndarray
- g_i
1-dimensional array of size N (=no of windows) containing WHAM-computed free energies for each window.
- Type
ndarray
Caution
All data must be at the same temperature.
Note
Binless WHAM handles empty bins when computing free energy profiles by setting bin free energy to inf.
Example
>>> calc = Calc1D() >>> status = calc.compute_point_weights(x_l, ...) >>> status True >>> G_l = calc.G_l >>> g_i = calc.g_i >>> betaF_x, _ = calc.bin_betaF_profile(x_l, G_l, x_bin, ...) >>> betaF_xy, _ = calc.bin_2D_betaF_profile(x_l, y_l, G_l, x_bin, y_bin, ...)
For comprehensive examples, check out the Jupyter notebooks in the examples/ folder.
- bin_2D_betaF_profile()
Bins weights corresponding to each sample into a 2D free energy profile in order parameters x_l (which is biased) and y_l (a related unbiased order parameter). If point weights G_l are not passed as an argument, then the computed WHAM weights are used for binning. You can pass custom weights G_l to compute reweighted free energy profiles.
Caution
This calculation uses the order parameter samples [self.x_l]. These will be available if you have called the compute function compute_point_weights or the main API call compute_betaF_profile. If you haven’t done so, you must initialize the Calc1D object’s x_l variable before calling this function.
- Parameters
y_l (ndarray) – Second dimension order parameter values.
x_bin (list) – Array of x-bin left-edges/centers of length M_x.
y_bin (list) – Array of y-bin left-edges/centers of length M_y.
G_l (ndarray) – Array of weights corresponding to each data point/order parameter x_l.
x_bin_style (string) – ‘left’ or ‘center’.
y_bin_style (string) – ‘left’ or ‘center’.
- Returns
- tuple(betaF_bin, tuple(betaF_bin_counts, delta_x_bin, delta_y_bin))
betaF_bin (ndarray): 2-D free energy profile of shape (M_x, M_y), binned using x_bin (1st dim) and y-bin (2nd dim).
betaF_bin_counts (ndarray): 2-D bin counts of shape (M_x, M_y)
delta_x_bin: Array of length M_x containing bin intervals along x.
delta_y_bin: Array of length M_y containing bin intervals along y.
- bin_betaF_profile()
Bins weights corresponding to each sample into a 1D free energy profile. If point weights G_l are not passed as an argument, then the computed WHAM weights are used for binning. You can pass custom weights G_l to compute reweighted free energy profiles.
Caution
This calculation uses the order parameter samples [self.x_l]. These will be available if you have called the compute function compute_point_weights or the main API call compute_betaF_profile. If you haven’t done so, you must initialize the Calc1D object’s x_l variable before calling this function.
- Parameters
x_bin (list) – Array of bin left-edges/centers of length M. Used only for computing final PMF.
G_l (ndarray) – Array of weights corresponding to each data point/order parameter x_l.
bin_style (string) – ‘left’ or ‘center’.
- Returns
Free energy profile, binned as per x_bin. betaF_bin_counts (ndarray): Bin counts
- Return type
betaF_bin (ndarray)
- bin_second_betaF_profile()
Bins weights corresponding to each sample into a into a 2D free energy profile in x_l (which is biased) and y_l (a related unbiased order parameter), then integrates out x_l to get a free energy profile in y_l. If point weights G_l are not passed as an argument, then the computed WHAM weights are used for binning. You can pass custom weights G_l to compute reweighted free energy profiles.
Caution
This calculation uses the order parameter samples [self.x_l]. These will be available if you have called the compute function compute_point_weights or the main API call compute_betaF_profile. If you haven’t done so, you must initialize the Calc1D object’s x_l variable before calling this function.
- Parameters
y_l (ndarray) – Second dimension order parameter values.
x_bin (list) – Array of x-bin left-edges/centers of length M_x.
y_bin (list) – Array of y-bin left-edges/centers of length M_y.
G_l (ndarray) – Array of weights corresponding to each data point/order parameter x_l.
x_bin_style (string) – ‘left’ or ‘center’.
y_bin_style (string) – ‘left’ or ‘center’.
- Returns
- Free energy profile of length M_y,
binned as per y-bin (2nd dim).
- Return type
betaF_y (ndarray)
- compute_betaF_profile()
Computes the binned free energy profile and window total free energies.
- Parameters
x_it (list) – Nested list of length S, x_it[i] is an array containing timeseries data from the i’th window.
x_bin (list) – Array of bin left-edges/centers of length M. Used only for computing final PMF.
u_i (list) – List of length S, u_i[i] is the umbrella potential function u_i(x) acting on the i’th window.
beta – beta, in inverse units to the units of u_i(x).
bin_style (string) – ‘left’ or ‘center’.
solver (string) – Solution technique to use [‘log-likelihood’, ‘self-consistent’, default=’log-likelihood’].
**solverkwargs – Arguments for solver.
- Returns
Free energy profile, binned as per x_bin. betaF_bin_counts (ndarray): Number of points in each bin. status (bool): Solver status.
- Return type
betaF_bin (ndarray)
- class WHAM.binless.CalcBase
Bases:
objectClass containing basic binless WHAM solver functionality.
- x_l
1-dimensional or 2-dimensional array containing unrolled order parameter which is being biased.
- Type
ndarray
- G_l
1-dimensional array containing WHAM-computed weights corresponding to each (unrolled) order parameter.
- Type
ndarray
- g_i
1-dimensional array of size N (=no of windows) containing WHAM-computed free energies for each window.
- Type
ndarray
Caution
All data must be at the same temperature.
Note
Binless WHAM handles empty bins when computing free energy profiles by setting bin free energy to inf.
- NLL()
Computes the negative log-likelihood objective function to minimize.
\(A = \frac{1}{N_{tot}} \left( \sum_{l=1}^{N_{tot}} \log \sum_{i=1}^{S} \frac{N_i}{N_{tot}} e^{g_i + W_{il}} - \sum_{i=1}^{S} N_i g_i \right)\)
- Parameters
g_i (ndarray of shape (S,)) – 0, 1, 2, …, S-1.
N_i (ndarray of shape (S,)) – Array of total sample counts for the windows 0, 1, 2, …, S-1.
W_il (ndarray of shape (S, Ntot)) – Array of weights, W_il = -beta * U_i(x_l) for the windows 0, 1, 2, …, S-1.
autograd (bool) – Use autograd to compute gradient.
- Returns
Negative log-likelihood objective function.
- Return type
A(g) (np.float)
- check_data()
Verifies that x_l is not None, else raises RuntimeError.
- Raises
RuntimeError –
- check_weights()
Verifies that g_i and G_l are not None, else raises RuntimeError.
- Raises
RuntimeError –
- compute_point_weights()
Computes WHAM weights corresponding to each order parameter sample and total window free energies. This is the main computation call for the Calc1D class.
- Parameters
x_l (ndarray) – Array containing each sample (unrolled).
N_i (ndarray) – Counts for each window.
u_i (list) – List of length S, u_i[i] is the umbrella potential function u_i(x) acting on the i’th window.
beta – beta, in inverse units to the units of u_i(x).
solver (string) – Solution technique to use [‘log-likelihood’, ‘self-consistent’, default=’log-likelihood’].
**solverkwargs – Arguments for solver.
- grad_NLL()
Computes the gradient of the negative log likelihood objective function w.r.t g_i.
\(\frac{\partial A}{\partial g_i} = \frac{1}{N_{tot}} \left( \sum_{l=1}^{N_{tot}} \frac{ N_i e^{g_i + W_{il}} }{ \sum_{i=1}^{S} N_i e^{g_i + W_{il}} } - N_i \right)\)
To use the log-sum-exp trick, this is computed as:
\(\frac{\partial A}{\partial g_i} = \frac{1}{N_{tot}} \left( e^{\Theta_i} - N_i \right)\)
\(\Theta_i = \log \sum_{l=1}^{N_{tot}} N_i e^{g_i + W_{il} - \Omega_l}\)
\(\Omega_l = \log \sum_{i=1}^{S} N_i e^{g_i + W_{il}}\)
- Parameters
g_i (ndarray of shape (S,)) – 0, 1, 2, …, S-1.
N_i (ndarray of shape (S,)) – Array of total sample counts for the windows 0, 1, 2, …, S-1.
W_il (ndarray of shape (S, Ntot)) – Array of weights, W_il = -beta * U_i(x_l) for the windows 0, 1, 2, …, S-1.
autograd (bool) – Use autograd to compute gradient (False). (Note: This is not evaluated.)
- Returns
Gradient of negative log-likelihood objective function w.r.t g_i.
- Return type
grad_A(g) (np.float)
- minimize_NLL_solver()
Computes optimal g_i by minimizing the negative log-likelihood for jointly observing the bin counts in the independent windows in the dataset.
Note
Any optimization method which scipy.optimize supports can be used. L-BFGS-B is used by default. Gradient information required for L-BFGS-B can either be computed analytically or using autograd.
- Parameters
N_i (ndarray of shape (S,)) – Array of total sample counts for the windows 0, 1, 2, …, S-1.
W_il (ndarray of shape (S, Ntot)) – Array of weights, W_il = -beta * U_i(x_l) for the windows 0, 1, 2, …, S-1.
g_i (ndarray of shape (S,)) – Total free energy initial guess.
opt_method (string) – Optimization algorithm to use (default: L-BFGS-B).
logevery (int) – Interval to log negative log-likelihood (default=100000000, i.e. no logging).
autograd (bool) – Use autograd to compute gradient.
- Returns
- tuple(bG_l, g_i, status)
bG_l (ndarray of shape (Ntot,)): Free energy for each sample point,
g_i (ndarray of shape (S,)): Total free energy for each window,
status (bool): Solution status.
- reweight()
Reweights sample weights to a biased ensemble. Does not change computed WHAM weights.
Caution
This is a post-processing calculation, and needs to be performed after computing weights through compute_point_weights or through the main API call compute_betaF_profile.
- Parameters
beta – beta: beta, in inverse units to the units of u_i(x).
u_bias – Biasing function applied to order parameter x_l.
- Returns
Biased weights for each data point x_l.
- Return type
G_l_bias (ndarray)
- self_consistent_solver()
Computes optimal parameters g_i by solving the coupled MBAR equations self-consistently until convergence.
- Parameters
N_i (ndarray of shape (S,)) – Array of total sample counts for the windows 0, 1, 2, …, S-1.
W_il (ndarray of shape (S, Ntot)) – Array of weights, W_il = -beta * U_i(x_l) for the windows 0, 1, 2, …, S-1.
g_i (ndarray of shape (S,)) – Total free energy initial guess.
tol (float) – Relative tolerance to stop solver iterations at (defaul=1e-7).
maxiter (int) – Maximum number of iterations to run solver for (default=100000).
logevery (int) – Interval to log self-consistent solver error (default=100000000, i.e. no logging).
- Returns
- tuple(bG_l, g_i, status)
bG_l (ndarray of shape (Ntot,)): Free energy for each sample point,
g_i (ndarray of shape (S,)): Total free energy for each window,
status (bool): Solution status.
- class WHAM.binless.CalcDD
Bases:
WHAM.binless.CalcBaseClass containing methods to compute free energy profiles from D-dimensional umbrella sampling data (i.e. data from biasing D order parameters) using binless WHAM.
- x_l
NxD matrix containing unrolled order parameter which is being biased, where (D=no of dimensions, N=no of data points).
- Type
ndarray
- G_l
1-dimensional array of length N (=no of data points) containing WHAM-computed weights corresponding to each (unrolled) order parameter.
- Type
ndarray
- g_i
1-dimensional array of size S (=no of windows) containing WHAM-computed free energies for each window.
- Type
ndarray
Caution
All data must be at the same temperature.
Note
Binless WHAM handles empty bins when computing free energy profiles by setting bin free energy to inf.
Example
>>> calc = CalcDD() >>> status = calc.compute_point_weights(X_l, ...) >>> status True >>> G_l = calc.G_l >>> g_i = calc.g_i >>> betaF_x, _ = calc.bin_betaF_profile(x_l, G_l, X_bin, ...)
For comprehensive examples, check out the Jupyter notebooks in the examples/ folder.
- betaF_stat()
- bin_betaF_profile()
Bins weights corresponding to each sample into a D-dimensional free energy profile. If point weights G_l are not passed as an argument, then the computed WHAM weights are used for binning. You can pass custom weights G_l to compute reweighted free energy profiles.
Caution
This calculation uses the order parameter samples [self.x_l]. These will be available if you have called the compute function compute_point_weights or the main API call compute_betaF_profile. If you haven’t done so, you must initialize the Calc1D object’s x_l variable before calling this function.
Important
The bins in each dimension must be uniform (unlike in Calc1D).
- Parameters
x_bin_list (list of ndarrays of dimensions (Nbin,)) – List of arrays of bin left edges/bin centers of length M. Used only for computing final PMF.
G_l (ndarray) – Array of weights corresponding to each data point/order parameter x_l.
in_style (string) – ‘left’ or ‘center’.
- Returns
Free energy profile, binned as per x_bin. betaF_bin_counts (ndarray): Bin counts
- Return type
betaF_bin (ndarray)
- compute_betaF_profile()
Computes the binned free energy profile and window total free energies.
- Parameters
x_it (list) – Nested list of length S, x_it[i] is an array containing timeseries data from the i’th window.
x_bin_list (list) – List of arrays of bin left edges/bin centers of length M. Used only for computing final PMF.
u_i (list) – List of length S, u_i[i] is the umbrella potential function u_i(x) acting on the i’th window.
beta – beta, in inverse units to the units of u_i(x).
bin_style (string) – ‘left’ or ‘center’.
solver (string) – Solution technique to use [‘log-likelihood’, ‘self-consistent’, default=’log-likelihood’].
**solverkwargs – Arguments for solver.
- Returns
Free energy profile, binned as per x_bin. betaF_bin_counts (ndarray): Number of points in each bin. status (bool): Solver status.
- Return type
betaF_bin (ndarray)
WHAM.binned module
- Implementation of binned WHAM using
Negative log-likelihood maximization as described in Zhu, F., & Hummer, G. (2012) with automatic differentiation.
Self-consistent iteration.
- class WHAM.binned.Calc1D
Bases:
objectClass containing methods to compute free energy profiles from 1D umbrella sampling data (i.e. data from biasing a single order parameter) using binned WHAM.
- x_l
1-dimensional array of bin left-edges/centers
- Type
ndarray
- betaF_l
1-dimensional array of size M (=no of bins) containing computed consensus free energy profile
- Type
ndarray
- g_i
1-dimensional array of size N (=no of windows) containing free energies for each window
- Type
ndarray
Caution
All data must be at the same temperature.
Notes
- If you want to compute 2D free energy profiles in x_l and a second (related, unbiased) order parameter y_l,
use binless WHAM instead, or use 2D binned WHAM.
Binned WHAM fails if bins are empty.
While binless WHAM is generally more accurate, binned WHAM is generally faster.
Example
>>> calc = Calc1D() >>> status = calc.compute_betaF_profile(...) >>> status True >>> betaF_l = calc.betaF_l >>> g_i = calc.g_i
For comprehensive examples, check out the Jupyter notebooks in the examples/ folder.
- NLL()
Computes the negative log-likelihood objective function to minimize.
- Parameters
g_i (np.array of shape (S,)) – 0, 1, 2, …, S-1.
N_i (np.array of shape (S,)) – Array of total sample counts for the windows 0, 1, 2, …, S-1.
M_l (np.array of shape (M,)) – Array of total bin counts for each bin.
W_il (np.array of shape (S, M)) – Array of weights, W_il = -beta * U_i(x_l) for the windows 0, 1, 2, …, S-1.
autograd (bool) – Use autograd to compute gradient.
- Returns
Negative log-likelihood objective function.
- Return type
A(g) (np.float)
- compute_betaF_profile()
Computes the binned free energy profile and window total free energies. Raises an Exception if any of the bins are empty.
- Parameters
x_it (list) – Nested list of length S, x_it[i] is an array containing timeseries data from the i’th window.
x_l (np.array) – Array of bin left-edges/centers of length M.
u_i (list) – List of length S, u_i[i] is the umbrella potential function u_i(x) acting on the i’th window.
beta – beta, in inverse units to the units of u_i(x).
bin_style (string) – ‘left’ or ‘center’.
solver (string) – Solution technique to use [‘log-likelihood’, ‘self-consistent’, ‘debug’, default=’log-likelihood’].
scale_stat_ineff (boolean) – Compute and scale bin count by statistical inefficiency (default=False).
**solverkwargs – Arguments for solver
- Returns
Solver status.
- Return type
status (bool)
- minimize_NLL_solver()
Computes optimal g_i by minimizing the negative log-likelihood for jointly observing the bin counts in the independent windows in the dataset.
Note
Any optimization method supported by scipy.optimize can be used. L-BFGS-B is used by default. Gradient information required for L-BFGS-B can either be computed analytically or using autograd.
- Parameters
N_i (np.array of shape (S,)) – Array of total sample counts for the windows 0, 1, 2, …, S-1.
M_l (np.array of shape (M,)) – Array of total bin counts for each bin.
W_il (np.array of shape (S, M)) – Array of weights, W_il = -beta * U_i(x_l) for the windows 0, 1, 2, …, S-1.
g_i (np.array of shape (S,)) – Total free energy initial guess.
opt_method (string) – Optimization algorithm to use (default: L-BFGS-B).
logevery (int) – Interval to log negative log-likelihood (default=0, i.e. no logging).
autograd (bool) – Use autograd to compute gradient.
- Returns
- tuple(g_i, status)
g_i: np.array of shape(S,))
status (bool): Solution status.
- self_consistent_solver()
Computes optimal parameters g_i by solving the coupled WHAM equations self-consistently until convergence.
- Parameters
N_i (np.array of shape (S,)) – Array of total sample counts for the windows 0, 1, 2, …, S-1.
M_l (np.array of shape (M,)) – Array of total bin counts for each bin.
W_il (np.array of shape (S, M)) – Array of weights, W_il = -beta * U_i(x_l) for the windows 0, 1, 2, …, S-1.
g_i (np.array of shape (S,)) – Total free energy initial guess.
tol (float) – Relative tolerance to stop solver iterations at (defaul=1e-7).
maxiter (int) – Maximum number of iterations to run solver for (default=100000).
logevery (int) – Interval to log self-consistent solver error (default=0, i.e. no logging).
- Returns
- tuple(g_i, status)
g_i (np.array of shape(S,)),
status (bool): Solution status.
- class WHAM.binned.CalcDD
Bases:
objectClass containing methods to compute free energy profiles from D-dimensional umbrella sampling data (i.e. data from biasing a D order parameters) using binless WHAM.
- X_l
DxM matrix containing unrolled order parameter which is being biased, where (D=no of dimensions, M=no of data points).
- Type
ndarray
- G_l
1-dimensional array of length M (=no of data points) containing WHAM-computed weights corresponding to each (unrolled) order parameter.
- Type
ndarray
- g_i
1-dimensional array of size N (=no of windows) containing WHAM-computed free energies for each window.
- Type
ndarray
Caution
All data must be at the same temperature.
Note
Binless WHAM handles empty bins when computing free energy profiles by setting bin free energy to inf.
Example
>>> calc = CalcND() >>> status = calc.compute_point_weights(X_l, ...) >>> status True >>> G_l = calc.G_l >>> g_i = calc.g_i >>> betaF_x, _ = calc.bin_betaF_profile(X_l, G_l, X_bin, ...)
For comprehensive examples, check out the Jupyter notebooks in the examples/ folder.
- X_l
- compute_betaF_profile()
WHAM.statistics module
Functions for checking consistency of WHAM calculations.
- WHAM.statistics.D_KL(betaF_P, betaF_Q, delta_x_bin)
Computes the KL divergence between two probability distributions P and Q corresponding to free energy profiles betaF_P and betaF_Q.
- Parameters
betaF_P (ndarray) – Free energy profile of length M.
betaF_Q (ndarray) – Free energy profile of length M.
delta_x_bin (ndarray) – Bin interval, length M.
- Returns
KL divergence (float)
- WHAM.statistics.D_KL_DD(betaF_P, betaF_Q, delta_X_bin)
- WHAM.statistics.binless_KLD_reweighted_win_betaF(calc, x_it, x_bin, u_i, beta, bin_style='left')
Computes the KL divergence between probability distributions corresponding to (sampled) bias free energy profile and reweighted biased free energy profile (constructed by reweighting binless WHAM profile) for each window.
- Parameters
calc (WHAM.binless.Calc1D) – Binless Calc1D object, with weights pre-computed.
x_it (list) – Nested list of length S, x_it[i] is an array containing timeseries data from the i’th window.
x_bin (list) – Array of bin left-edges/centers of length M. Used only for computing final PMF.
u_i (list) – List of length S, u_i[i] is the umbrella potential function u_i(x) acting on the i’th window.
beta – beta, in inverse units to the units of u_i(x).
bin_style (string) – ‘left’ or ‘center’.
- Returns
Array of length S containing KL divergence for each window
- Return type
D_KL_i (ndarray)
- WHAM.statistics.binless_KLD_reweighted_win_betaF_DD(calc, X_it, X_bin, u_i, beta, bin_style='left')
- WHAM.statistics.binless_reweight_phi_ensemble(calc, phi_vals, beta)
Computes <x_l> and <dx_l^2> v/s phi by reweighting the underlying free energy profile in calc to the phi-ensemble specified by phi_vals and beta.
- Parameters
calc (WHAM.binless.Calc1D) – Binless Calc1D object, with weights pre-computed.
phi_vals (ndarray) – Array of phi values to compute <x_l> at, units of phi (e.g. kJ/mol).
beta – beta, in inverse units to the units of u_i(x).
bin_style (string) – ‘left’ or ‘center’.
- Returns
Array containing <x_l> corresponding to each value of phi. x_var_vals (ndarray): Array containing <dx_l^2> corresponding to each value of phi.
- Return type
x_avg_vals (ndarray)
- WHAM.statistics.binless_reweighted_win_betaF(calc, x_bin, u_i, beta, bin_style='left')
Computes free energy profiles for each window i by reweighting the entire free energy profile to the biased ensemble with bias potential specified by u_i[i].
- Parameters
calc (WHAM.binless.Calc1D) – Binless Calc1D object, with weights pre-computed.
x_bin (list) – Array of bin left-edges/centers of length M. Used only for computing final PMF.
u_i (list) – List of length S, u_i[i] is the umbrella potential function u_i(x) acting on the i’th window.
beta – beta, in inverse units to the units of u_i(x).
bin_style (string) – ‘left’ or ‘center’.
- Returns
2-D array of biased free energies for each window, of shape (S, M)
- Return type
betaF_il_reweight (np.array)
- WHAM.statistics.binless_reweighted_win_betaF_DD(calc, X_bin, u_i, beta, bin_style='left')
- WHAM.statistics.binned_KLD_reweighted_win_betaF(x_it, x_bin, betaF_bin, u_i, beta, bin_style='left')
Computes the KL divergence between probability distributions corresponding to (sampled) bias free energy profile and reweighted biased free energy profile (constructed by reweighting binned WHAM profile) for each window.
- Parameters
x_it (list) – Nested list of length S, x_it[i] is an array containing timeseries data from the i’th window.
x_bin (list) – Array of bin left-edges/centers of length M.
betaF_bin (list) – Array of (consensus) free energies of length M.
u_i (list) – List of length S, u_i[i] is the umbrella potential function u_i(x) acting on the i’th window.
beta – beta, in inverse units to the units of u_i(x).
bin_style (string) – ‘left’ or ‘center’.
Caution
For the best results, use the same bin style as the consensus free energy profile.
- Returns
Array of length S containing KL divergence for each window
- Return type
D_KL_i (ndarray)
- WHAM.statistics.binned_KLD_reweighted_win_betaF_DD(X_it, X_bin, betaF_bin, u_i, beta, bin_style='left')
- WHAM.statistics.binned_reweight_phi_ensemble(x_bin, betaF_bin, phi_vals, beta)
Computes <x_l> and <dx_l^2> v/s phi by reweighting the free energy profile betaF_bin to the phi-ensemble specified by phi_vals and beta.
- Parameters
x_bin (list) – Array of bin left-edges/centers of length M.
phi_vals (ndarray) – Array of phi values to compute <x_l> at, units of phi (e.g. kJ/mol).
betaF_bin (list) – Array of free energies of length M.
beta – beta, in inverse units to the units of u_i(x).
- Returns
Array containing <x_l> corresponding to each value of phi. x_var_vals (ndarray): Array containing <dx_l^2> corresponding to each value of phi.
- Return type
x_avg_vals (ndarray)
- WHAM.statistics.binned_reweighted_win_betaF(x_bin, betaF_bin, u_i, beta)
Computes free energy profile for each window i by reweighting the entire free energy profile betaF_bin at x_bin to the biased ensemble with bias potential specified by u_i[i].
Note
This function can be also be used to reweight free energy profiles constructed using binless WHAM.
- Parameters
x_bin (list) – Array of bin left-edges/centers of length M.
betaF_bin (list) – Array of free energies of length M.
u_i (list) – List of length S, u_i[i] is the umbrella potential function u_i(x) acting on the i’th window.
beta – beta, in inverse units to the units of u_i(x).
- Returns
2-D array of biased free energies for each window, of shape (S, M)
- Return type
betaF_il_reweight (np.array)
- WHAM.statistics.binned_reweighted_win_betaF_DD(X_bin, betaF_bin, u_i, beta)
- WHAM.statistics.win_betaF(x_it, x_bin, u_i, beta, bin_style='left', scale_stat_ineff=False)
Computes biased free energy profiles for each umbrella window from raw data.
- Parameters
x_it (list) – Nested list of length S, x_it[i] is an array containing timeseries data from the i’th window.
x_bin (list) – Array of bin left-edges/centers of length M. Used only for computing final PMF.
u_i (list) – List of length S, u_i[i] is the umbrella potential function u_i(x) acting on the i’th window.
beta – beta, in inverse units to the units of u_i(x).
bin_style (string) – ‘left’ or ‘center’.
- Returns
- tuple(betaF_il, delta_x_bin)
betaF_il (np.array): 2-D array of biased free energies for each window, of shape (S, M)
delta_x_bin (np.array): Bin interval
- WHAM.statistics.win_betaF_DD(X_it, X_bin, u_i, beta, bin_style='left', scale_stat_ineff=False)