API#
- class bayesmbar.BayesBAR(energy: ndarray, num_conf: ndarray, sample_size: int = 1000, method: str = 'Newton', verbose=False)[source]#
Bases:
object
Bayesian Bennett acceptance ratio method
Initialize the BayesBAR class.
- Parameters:
energy (ndarray) – An energy matrix in reduced units. Its size should be 2xN, where N is the total number of samples from the two states.
num_conf (ndarray) – Number of configurations in each state. Its size should be (2,).
sample_size (int, optional) – The number of samples from the posterior distribution. Defaults to 1000.
method (str, optional) – Optimization method for finding the mode. Options are “Newton” or “L-BFGS-B”. Defaults to “Newton”.
verbose (bool, optional) – Whether to print running information. Defaults to False.
- property DeltaF_mean: ndarray#
The posterior mean of the free energy difference.
- property DeltaF_mode: ndarray#
The posterior mode of the free energy difference.
- property DeltaF_samples: ndarray#
The samples from the posterior distribution of the free energy difference.
- property DeltaF_std: ndarray#
The posterior standard deviation of the free energy difference.
- class bayesmbar.BayesMBAR(energy: ndarray, num_conf: ndarray, prior: Literal['uniform', 'normal'] = 'uniform', mean: Literal['constant', 'linear', 'quadratic'] = 'constant', kernel: Literal['SE', 'Matern52', 'Matern32', 'RQ'] = 'SE', state_cv: ndarray = None, sample_size: int = 1000, warmup_steps: int = 500, optimize_steps: int = 10000, verbose: bool = True, random_seed: int = 0, method: Literal['Newton', 'L-BFGS-B'] = 'Newton')[source]#
Bases:
object
Bayesian Multistate Bennett Acceptance Ratio (BayesMBAR) method
- Parameters:
energy (np.ndarray) – Energy matrix of shape (m, n), where m is the number of states and n is the number of configurations.
num_conf (np.ndarray) – Number of configurations in each state. It is a 1D array of length m.
prior (str, optional) – Prior distribution of dF. It can be either “uniform” or “normal”. Defaults to “uniform”.
mean (str, optional) – Mean function of the prior. It can be either “constant”, “linear”, or “quadratic”. Defaults to “constant”.
kernel (str, optional) – Kernel function of the prior. It can be either “SE”, “Matern52”, “Matern32”, or “RQ”. Defaults to “SE”.
state_cv (np.ndarray, optional) – State collective variables. It is a 2D array of shape (n, d), where n is the number of configurations and d is the dimension of the collective variables. Defaults to None.
sample_size (int, optional) – Number of samples drawn from the posterior distribution. Defaults to 1000.
warmup_steps (int, optional) – Number of warmup steps used to find the step size and mass matrix of the NUTS sampler. Defaults to 500.
optimize_steps (int, optional) – Number of optimization steps used to learn the hyperparameters when normal priors are used. Defaults to 10000.
verbose (bool, optional) – Whether to print the progress bar for the optimization and sampling. Defaults to True.
random_seed (int, optional) – Random seed. Defaults to 0.
- property DeltaF_mean: ndarray[tuple[int, ...], dtype[_ScalarType_co]]#
The posterior mean of free energy difference between states. DeltaF_mean[i,j] is the free energy difference between state \(j\) and state \(i\), i.e., DeltaF_mean[i,j] = F_mean[j] - F_mean[i].
- property DeltaF_mode: ndarray[tuple[int, ...], dtype[_ScalarType_co]]#
The posterior mode estimate of free energy difference between states. DeltaF_mode[i,j] is the free energy difference between state \(j\) and state \(i\), i.e., DeltaF_mode[i,j] = F_mode[j] - F_mode[i].
- property DeltaF_std: ndarray[tuple[int, ...], dtype[_ScalarType_co]]#
The posterior standard deviation of free energy difference between states. DeltaF_std[i,j] is the posterior standard deviation of the free energy difference between state \(j\) and state \(i\),
- property F_cov: ndarray[tuple[int, ...], dtype[_ScalarType_co]]#
The posterior covariance matrix of the free energies of the states under the constraints that \(\\sum_{k=1}^{M} N_k * F_k = 0\), where \(N_k\) and \(F_k\) are the number of conformations and the free energy of the k-th state, respectively.
- property F_mean: ndarray[tuple[int, ...], dtype[_ScalarType_co]]#
The posterior mean of the free energies of the states under the constraints that \(\\sum_{k=1}^{M} N_k * F_k = 0\), where \(N_k\) and \(F_k\) are the number of conformations and the free energy of the k-th state, respectively.
- property F_mode: ndarray[tuple[int, ...], dtype[_ScalarType_co]]#
The posterior mode estimate of the free energies of the states under the constraints that \(\sum_{k=1}^{M} N_k * F_k = 0\), where \(N_k\) and \(F_k\) are the number of conformations and the free energy of the k-th state, respectively.
- property F_samples: ndarray[tuple[int, ...], dtype[_ScalarType_co]]#
The samples of the free energies of the states from the posterior distribution under the constraints that \(\sum_{k=1}^{M} N_k * F_k = 0\), where \(N_k\) and \(F_k\) are the number of conformations and the free energy of the k-th state, respectively.
- property F_std: ndarray[tuple[int, ...], dtype[_ScalarType_co]]#
The posterior standard deviation of the free energies of the states under the constraints that \(\\sum_{k=1}^{M} N_k * F_k = 0\), where \(N_k\) and \(F_k\) are the number of conformations and the free energy of the k-th state, respectively.
- class bayesmbar.CBayesMBAR(energies: List[ndarray], nums_conf: List[ndarray], identical_states: List[List[Tuple[int, int]]], sample_size: int = 1000, warmup_steps: int = 500, method: str = 'Newton', random_seed: int = None, verbose: bool = True)[source]#
Bases:
object
Coupled BayesMBAR
- Parameters:
energies (List[ndarray]) – A list of energies for all coupled MBAR systems. The energies should be in the unit of kT.
nums_conf (List[ndarray]) – A list of the number of configurations for all coupled MBAR systems.
identical_states (List[List[(int, int)]]) – A list of identical states. Each element in the outer list is a list of tuples and represents a group of states that are identical to each other. States are represented by a tuple where the first element is the index of a MBAR system and the second element is the index of the state in that MBAR system. For example, identical_states = [[(0, 3), (1, 0)], [(1, 3), (2, 0), (3, 0)]] means that state 3 in system 0 is identical to state 0 in system 1, and state 3 in system 1 is identical to state 0 in system 2, and state 0 in system 3.
sample_size (int, optional) – The number of samples to draw from the likelihood. The default is 1000.
warmup_steps (int, optional) – The number of warmup steps for the HMC sampler. The default is 500.
method (str, optional) – The optimization method for finding the mode of the likelihood. Options are “Newton” or “L-BFGS-B”. The default is “Newton”.
random_seed (int, optional) – The random seed. The default is None, which means the random seed is generated from the current time.
verbose (bool, optional) – Whether to print out the progress of the sampling. The default is True.
- property DeltaF_mean: List[ndarray]#
The mean of free energy differences between all pairs of states in every MBAR system.
- property DeltaF_mode: List[ndarray]#
The mode of free energy differences between all pairs of states in every MBAR system.
- property DeltaF_std: List[ndarray]#
The standard deviation of free energy differences between all pairs of states in every MBAR system.
- property F_mean: List[ndarray]#
The mean of free energies of all states in all MBAR systems. The free energy of state 0 in each system is set to 0.
- property F_mode: List[ndarray]#
The mode of free energies of all states in all MBAR systems. The free energy of state 0 in each system is set to 0.
- property F_samples: List[ndarray]#
The samples of free energies of all states in all MBAR systems. The free energy of state 0 in each system is set to 0.
- class bayesmbar.FastMBAR(energy: ndarray[tuple[int, ...], dtype[float64]], num_conf: ndarray[tuple[int, ...], dtype[int64]], bootstrap: bool = False, bootstrap_block_size: int = 3, bootstrap_num_rep: int = 100, verbose: bool = False, method: str = 'Newton')[source]#
Bases:
object
The FastMBAR class is initialized with an energy matrix and an array of num of conformations. The corresponding MBAR equation is solved in the constructor. Therefore, the relative free energies of states used in the energy matrix is calculated in the constructor. The method calculate_free_energies_for_perturbed_states can be used to calculated the relative free energies of perturbed states.
Initializer for the class FastMBAR
- Parameters:
energy (2D ndarray) – It has a size of (M x N), where M is the number of states and N is the total number of conformations. The entry energy[i,j] is the reduced (unitless) energy of conformation j in state i. If bootstrapping is used to calculate the uncertainty, the order of conformations matters. Conformations sampled from one state need to occpy a continuous chunk of collumns. Conformations sampled from state k need to occupy collumns to the left of conformations sampled from state l if k < l. If bootstrapping is not used, then the order of conformation does not matter.
num_conf (1D int ndarray) – It should have a size of M, where num_conf[i] is the num of conformations sampled from state i. Therefore, np.sum(num_conf) has to be equal to N. All entries in num_conf have to be strictly greater than 0.
bootstrap (bool, optional (default=False)) – If bootstrap is True, the uncertainty of the calculated free energies will be estimate using block bootstraping.
bootstrap_block_size (int, optional (default=3)) – block size used in block bootstrapping
bootstrap_num_rep (int, optional (default=100)) – number of repreats in block bootstrapping
verbose (bool, optional (default=False)) – if verbose is true, the detailed information of solving MBAR equations is printed.
method (str, optional (default="Newton")) – the method used to solve the MBAR equation. The default is Newton’s method.
- property DeltaF: ndarray#
Free energy difference between states. \(\mathrm{DeltaF}[i,j]\) is the free energy difference between state j and state i, i.e., \(\mathrm{DeltaF}[i,j] = F[j] - F[i]\) .
- property DeltaF_std: ndarray#
Standard deviation of the free energy difference between states. \(\mathrm{DeltaF_std}[i,j]\) is the standard deviation of the free energy difference \(\mathrm{DeltaF}[i,j]\).
- property F: ndarray#
Free energies of the states under the constraint \(\sum_{k=1}^{M} N_k * F_k = 0\), where \(N_k\) is the number of conformations sampled from state k.
- property F_cov: ndarray#
Covariance matrix of the free energies of the states under the constraint \(\sum_{k=1}^{M} N_k * F_k = 0\), where \(N_k\) is the number of conformations sampled from state k.
- property F_std: ndarray#
Standard deviation of the free energies of the states under the constraint \(\sum_{k=1}^{M} N_k * F_k = 0\), where \(N_k\) is the number of conformations sampled from state k.
- calculate_free_energies_of_perturbed_states(energy_perturbed: ndarray[tuple[int, ...], dtype[float64]]) dict [source]#
calculate free energies for perturbed states.
- Parameters:
energy_perturbed (2-D float ndarray with size of (L,N)) – each row of the energy_perturbed matrix represents a state and the value energy_perturbed[l,n] represents the reduced energy of the n’th conformation in the l’th perturbed state.
- Returns:
results – a dictionary containing the following keys:
F - the free energy of the perturbed states.
F_std - the standard deviation of the free energy of the perturbed states.
F_cov - the covariance between the free energies of the perturbed states.
DeltaF - \(\mathrm{DeltaF}[k,l]\) is the free energy difference between state \(k\) and state \(l\), i.e., \(\mathrm{DeltaF}[k,l] = F[l] - F[k]\) .
DeltaF_std - the standard deviation of the free energy difference.
- Return type:
dict
- property log_prob_mix: ndarray#
the log probability density of conformations in the mixture distribution.