FastMBAR API#

class FastMBAR.FastMBAR(energy: ndarray | Tensor, num_conf: ndarray | Tensor, cuda: bool = False, cuda_batch_mode: bool | None = None, bootstrap: bool = False, bootstrap_block_size: int = 3, bootstrap_num_rep: int = 100, verbose: bool = False, method: str = 'Newton')[source]#

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.

__init__(energy: ndarray | Tensor, num_conf: ndarray | Tensor, cuda: bool = False, cuda_batch_mode: bool | None = None, bootstrap: bool = False, bootstrap_block_size: int = 3, bootstrap_num_rep: int = 100, verbose: bool = False, method: str = 'Newton') None[source]#

Initializer for the class FastMBAR

Parameters:
  • energy (2D ndarray or 2D tensor) – 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 or 1D tensor) – 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.

  • cuda (bool, optional (default=False)) – If it is set to be True, then the calculation will be run on a graphical processing unit (GPU) using CUDA.

  • cuda_batch_mode (bool, optional (default=None)) – The batch mode is turned on when the size of the energy matrix is too large for the memory of a GPU. When cuda_batch_mode is True, the energy matrix will be split into multiple batches which are used sequentially. If cuda_batch_mode = None, it will be set automatically based on the size of energy and the memory of the avaible GPU device.

  • 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 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_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.

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 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 log_prob_mix: ndarray#

the log probability density of conformations in the mixture distribution.

calculate_free_energies_of_perturbed_states(energy_perturbed)[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