umbrella sampling module
- class umbrella_sampling.BaseUmbrellaSampling(file_dir, system, clean_build=False)
- wham_run(wham_dir, xmin, xmax, umbrella_stiff, n_bins, tol, n_boot, all_observables=False)
Run the weighted histogram analysis technique (Grossfield, Alan http://membrane.urmc.rochester.edu/?page_id=126)
- Parameters:
wham_dir (str): directory in which wham is complied using the install_wham.sh script.
xmin (str): smallest value of the order parameter to be sampled.
xmax (str): largest value of the order parameter to be sample.
umbrella_stiff (str): stiffness of the umbrella potential.
n_bins (str): number of bins of the resultant free energy profile.
tol (str): the tolerance of convergence for the WHAM technique.
n_boot (str): number of monte carlo bootstrapping steps to take.
- class umbrella_sampling.ComUmbrellaSampling(file_dir, system, clean_build=False)
- build_equlibration_runs(simulation_manager, n_windows, com_list, ref_list, stiff, xmin, xmax, input_parameters, observable=False, sequence_dependant=False, print_every=10000.0, name='com_distance.txt', continue_run=False, protein=None, force_file=None)
Build the equlibration run
- Parameters:
simulation_manager (SimulationManager): pass a simulation manager.
n_windows (int): number of simulation windows.
com_list (str): comma speperated list of nucleotide indexes.
ref_list (str): comma speperated list of nucleotide indexex.
stiff (float): stiffness of the umbrella potential.
xmin (float): smallest value of the order parameter to be sampled.
xmax (float): largest value of the order parameter to be sampled.
input_parameters (dict): dictonary of oxDNA parameters.
observable (bool): Boolean to determine if you want to print the observable for the equlibration simulations.
sequence_dependant (bool): Boolean to use sequence dependant parameters.
print_every (float): how often to print the order parameter.
name (str): depreciated.
continue_run (float): number of steps to continue umbrella sampling (i.e 1e7).
protein (bool): Use a protein par file and run with the ANM oxDNA interaction potentials.
force_file (bool): Add an external force to the umbrella simulations.
- build_production_runs(simulation_manager, n_windows, com_list, ref_list, stiff, xmin, xmax, input_parameters, observable=True, sequence_dependant=False, print_every=10000.0, name='com_distance.txt', continue_run=False, protein=None, force_file=None)
Build the production run
- Parameters:
simulation_manager (SimulationManager): pass a simulation manager.
n_windows (int): number of simulation windows.
com_list (str): comma speperated list of nucleotide indexes.
ref_list (str): comma speperated list of nucleotide indexex.
stiff (float): stiffness of the umbrella potential.
xmin (float): smallest value of the order parameter to be sampled.
xmax (float): largest value of the order parameter to be sampled.
input_parameters (dict): dictonary of oxDNA parameters.
observable (bool): Boolean to determine if you want to print the observable for the equlibration simulations, make sure to not have this as false for production run.
sequence_dependant (bool): Boolean to use sequence dependant parameters.
print_every (float): how often to print the order parameter.
name (str): depreciated.
continue_run (float): number of steps to continue umbrella sampling (i.e 1e7).
protein (bool): Use a protein par file and run with the ANM oxDNA interaction potentials.
force_file (bool): Add an external force to the umbrella simulations.
- umbrella_forces(com_list, ref_list, stiff, xmin, xmax, n_windows)
Build Umbrella potentials
- rate_umbrella_forces(com_list, ref_list, stiff, xmin, xmax, n_windows, starting_r0, steps)
Build Umbrella potentials
- class umbrella_sampling.WhamAnalysis(base_umbrella)
- run_wham(wham_dir, xmin, xmax, umbrella_stiff, n_bins, tol, n_boot)
Run Weighted Histogram Analysis Method on production windows.
- Parameters:
wham_dir (str): Path to wham executable. xmin (str): Minimum distance of center of mass order parameter in simulation units. xmax (str): Maximum distance of center of mass order parameter in simulation units. umbrella_stiff (str): The parameter used to modified the stiffness of the center of mass spring potential n_bins (str): number of histogram bins to use. tol (str): Convergence tolerance for the WHAM calculations. n_boot (str): Number of monte carlo bootstrapping error analysis iterations to preform.
- chunk_convergence_analysis(n_chunks)
Seperate your data into equal chunks
- data_truncated_convergence_analysis(data_added_per_iteration)
Seperate your data into equal chunks
- convergence_analysis(n_chunks, data_added_per_iteration, wham_dir, xmin, xmax, umbrella_stiff, n_bins, tol, n_boot)
Split your data into a set number of chunks to check convergence. If all chunks are the same your free energy profile if probabily converged Also create datasets with iteraterativly more data, to check convergence progress