simulation.rep_exch submodule
- cg_openmm.simulation.rep_exch.extract_trajectory(topology, reporter, state_index=None, replica_index=None, frame_begin=0, frame_stride=1, frame_end=-1)[source]
Internal function for extract trajectory (replica or state) from .nc file, Based on YANK extract_trajectory code.
- cg_openmm.simulation.rep_exch.get_minimum_energy_ensemble(topology, replica_energies, replica_positions, ensemble_size=5, file_name=None)[source]
Get an ensemble of low (potential) energy poses, and write the lowest energy structure to a PDB file if a file_name is provided.
- Parameters
topology (Topology()) – OpenMM Topology()
replica_energies (List( List( float * simtk.unit.energy for simulation_steps ) for num_replicas )) – List of dimension num_replicas X simulation_steps, which gives the energies for all replicas at all simulation steps
replica_positions (np.array( ( float * simtk.unit.positions for num_beads ) for simulation_steps )) – List of positions for all output frames for all replicas
file_name – Output destination for PDB coordinates of minimum energy pose, Default = None
- Returns
ensemble ( List() ) - A list of poses that are in the minimum energy ensemble.
- Example
>>> from foldamers.cg_model.cgmodel import CGModel >>> from cg_openmm.simulation.rep_exch import * >>> cgmodel = CGModel() >>> replica_energies,replica_positions,replica_state_indices = run_replica_exchange(cgmodel.topology,cgmodel.system,cgmodel.positions) >>> ensemble_size = 5 >>> file_name = "minimum.pdb" >>> minimum_energy_ensemble = get_minimum_energy_ensemble(cgmodel.topology,replica_energies,replica_positions,ensemble_size=ensemble_size,file_name=file_name)
- cg_openmm.simulation.rep_exch.make_replica_dcd_files(topology, timestep=Quantity(value=5, unit=femtosecond), time_interval=200, output_dir='output', output_data='output.nc', checkpoint_data='output_checkpoint.nc', frame_begin=0, frame_stride=1, center=False)[source]
Make dcd files from replica exchange simulation trajectory data.
- Parameters
topology – OpenMM Topology
timestep (Quantity() <http://docs.openmm.org/development/api-python/generated/simtk.unit.quantity.Quantity.html> float * simtk.unit) – Time step used in the simulation (default=5*unit.femtosecond)
time_interval (int) – frequency, in number of time steps, at which positions were recorded (default=200)
output_dir (str) – path to which we will write the output (default=’output’)
output_data (str) – name of output .nc data file (default=’output.nc’)
checkpoint_data (str) – name of checkpoint .nc data file (default=’output_checkpoint.nc’)
frame_begin (int) – Frame at which to start writing the dcd trajectory (default=0)
frame_stride (int) – advance by this many time intervals when writing dcd trajectories (default=1)
center (Boolean) – align all frames in the replica trajectories (default=False)
- cg_openmm.simulation.rep_exch.make_replica_pdb_files(topology, output_dir='output', output_data='output.nc', checkpoint_data='output_checkpoint.nc', frame_begin=0, frame_stride=1, center=False)[source]
Make pdb files from replica exchange simulation trajectory data.
- Parameters
topology – OpenMM Topology
output_dir (str) – path to which we will write the output (default=’output’)
output_data (str) – name of output .nc data file (default=’output.nc’)
checkpoint_data (str) – name of checkpoint .nc data file (default=’output_checkpoint.nc’)
frame_begin (int) – Frame at which to start writing the pdb trajectory (default=0)
frame_stride (int) – advance by this many frames when writing pdb trajectories (default=1)
center (Boolean) – align all frames in the replica trajectories (default=False)
- Returns
file_list ( List( str ) ) - A list of names for the files that were written
- cg_openmm.simulation.rep_exch.make_state_dcd_files(topology, timestep=Quantity(value=5, unit=femtosecond), time_interval=200, output_dir='output', output_data='output.nc', checkpoint_data='output_checkpoint.nc', frame_begin=0, frame_stride=1, center=True)[source]
Make dcd files by state from replica exchange simulation trajectory data. Note: these are discontinuous trajectories with constant temperature state.
- Parameters
topology – OpenMM Topology
timestep (Quantity() <http://docs.openmm.org/development/api-python/generated/simtk.unit.quantity.Quantity.html> float * simtk.unit) – Time step used in the simulation (default=5*unit.femtosecond)
time_interval (int) – frequency, in number of time steps, at which positions were recorded (default=200)
output_dir (str) – path to which we will write the output (default=’output’)
output_data (str) – name of output .nc data file (default=’output.nc’)
checkpoint_data (str) – name of checkpoint .nc data file (default=’output_checkpoint.nc’)
frame_begin (int) – Frame at which to start writing the dcd trajectory (default=0)
frame_stride (int) – advance by this many time intervals when writing dcd trajectories (default=1)
center (Boolean) – align the center of mass of each structure in the discontinuous state trajectory (default=True)
- cg_openmm.simulation.rep_exch.make_state_pdb_files(topology, output_dir='output', output_data='output.nc', checkpoint_data='output_checkpoint.nc', frame_begin=0, frame_stride=1, center=True)[source]
Make pdb files by state from replica exchange simulation trajectory data. Note: these are discontinuous trajectories with constant temperature state.
- Parameters
topology – OpenMM Topology
output_dir (str) – path to which we will write the output (default=’output’)
output_data (str) – name of output .nc data file (default=’output.nc’)
checkpoint_data (str) – name of checkpoint .nc data file (default=’output_checkpoint.nc’)
frame_begin (int) – Frame at which to start writing the pdb trajectory (default=0)
frame_stride (int) – advance by this many frames when writing pdb trajectories (default=1)
center (Boolean) – align the center of mass of each structure in the discontinuous state trajectory (default=True)
- Returns
file_list ( List( str ) ) - A list of names for the files that were written
- cg_openmm.simulation.rep_exch.plot_replica_exchange_energies(state_energies, temperature_list, series_per_page, time_interval=Quantity(value=1.0, unit=picosecond), time_shift=Quantity(value=0.0, unit=picosecond), file_name='rep_ex_ener.pdf', legend=True)[source]
Plot the potential energies for a batch of replica exchange trajectories
- Parameters
state_energies (List( List( float * simtk.unit.energy for simulation_steps ) for num_replicas )) – List of dimension num_replicas X simulation_steps, which gives the energies for all replicas at all simulation steps
temperature_list – List of temperatures for which to perform replica exchange simulations, default = [(300.0 * unit.kelvin).__add__(i * unit.kelvin) for i in range(-20,100,10)]
time_interval (SIMTK Unit()) – interval between energy exchanges.
time_shift – amount of time before production period to shift the time axis(default = 0)
file_name (str) – The pathname of the output file for plotting results, default = “replica_exchange_energies.png”
legend (Logical) – Controls whether a legend is added to the plot
- cg_openmm.simulation.rep_exch.plot_replica_exchange_energy_histograms(state_energies, temperature_list, file_name='rep_ex_ener_hist.pdf', legend=True)[source]
Plot the potential energies for a batch of replica exchange trajectories
- Parameters
state_energies (List( List( float * simtk.unit.energy for simulation_steps ) for num_replicas )) – List of dimension num_replicas X simulation_steps, which gives the energies for all replicas at all simulation steps
temperature_list – List of temperatures for which to perform replica exchange simulations, default = [(300.0 * unit.kelvin).__add__(i * unit.kelvin) for i in range(-20,100,10)]
file_name (str) – The pathname of the output file for plotting results, default = “replica_exchange_energies.png”
legend (Logical) – Controls whether a legend is added to the plot
- cg_openmm.simulation.rep_exch.plot_replica_exchange_summary(replica_states, temperature_list, series_per_page, time_interval=Quantity(value=1.0, unit=picosecond), time_shift=Quantity(value=0.0, unit=picosecond), file_name='rep_ex_states.pdf', legend=True)[source]
Plot the thermodynamic state assignments for individual temperature replicas as a function of the simulation time, in order to obtain a visual summary of the replica exchanges from a OpenMM simulation.
- Parameters
replica_states (List( List( float * simtk.unit.energy for simulation_steps ) for num_replicas )) – List of dimension num_replicas X simulation_steps, which gives the thermodynamic state indices for all replicas at all simulation steps
temperature_list – List of temperatures for which to perform replica exchange simulations, default = [(300.0 * unit.kelvin).__add__(i * unit.kelvin) for i in range(-20,100,10)]
time_interval – interval between energy exchanges.
time_shift – amount of time before production period to shift the time axis(default = 0)
file_name (str) – The pathname of the output file for plotting results, default = “replica_exchange_state_transitions.png”
legend (Logical) – Controls whether a legend is added to the plot
- cg_openmm.simulation.rep_exch.process_replica_exchange_data(output_data='output/output.nc', output_directory='output', series_per_page=4, write_data_file=True, plot_production_only=False, print_timing=False, equil_nskip=1, frame_begin=0, frame_end=-1)[source]
Read replica exchange simulation data, detect equilibrium and decorrelation time, and plot replica exchange results.
- Parameters
output_data (str) – path to output .nc file from replica exchange simulation, (default=’output/output.nc’)
output_directory (stry) – path to which output files will be written (default=’output’)
series_per_page (int) – number of replica data series to plot per pdf page (default=4)
write_data_file (Boolean) – Option to write a text data file containing the state_energies array (default=True)
plot_production_only (Boolean) – Option to plot only the production region, as determined from pymbar detectEquilibration (default=False)
equil_nskip (Boolean) – skip this number of frames to sparsify the energy timeseries for pymbar detectEquilibration (default=1) - this is used only when frame_begin=0 and the trajectory has less than 40000 frames.
frame_begin (int) – analyze starting from this frame, discarding all prior as equilibration period (default=0)
frame_end (int) – analyze up to this frame only, discarding the rest (default=-1).
- Returns
replica_energies ( Quantity() ( np.float( [number_replicas,number_simulation_steps] ), simtk.unit ) ) - The potential energies for all replicas at all (printed) time steps
replica_state_indices ( np.int64( [number_replicas,number_simulation_steps] ), simtk.unit ) - The thermodynamic state assignments for all replicas at all (printed) time steps
production_start ( int - The frame at which the production region begins for all replicas, as determined from pymbar detectEquilibration
sample_spacing ( int - The number of frames between uncorrelated state energies, estimated using heuristic algorithm )
n_transit ( np.float( [number_replicas] ) ) - Number of half-transitions between state 0 and n for each replica
mixing_stats ( tuple ( np.float( [number_replicas x number_replicas] ) , np.float( [ number_replicas ] ) , float( statistical inefficiency ) ) ) - transition matrix, corresponding eigenvalues, and statistical inefficiency
- cg_openmm.simulation.rep_exch.restart_replica_exchange(total_simulation_time, simulation_time_step=Quantity(value=5, unit=picosecond), exchange_frequency=200, output_data='output/output.nc')[source]
Restart an OpenMMTools replica exchange simulation using an OpenMM coarse grained model and output .nc files from the previous segment of the simulation.
- Parameters
total_simulation_time – Total run time for original + new simulation segments
simulation_time_step – Simulation integration time step (default=5*unit.picosecond)
exchange_frequency (int) – Number of time steps between replica exchange attempts (default=200)
output_data (str) – Path to the NETCDF file for previous segment of simulation - this will be appended to (default=”output/output.nc”)
- cg_openmm.simulation.rep_exch.run_replica_exchange(topology, system, positions, total_simulation_time=Quantity(value=1.0, unit=picosecond), simulation_time_step=None, temperature_list=None, friction=Quantity(value=1.0, unit=/picosecond), minimize=True, exchange_frequency=1000, output_data='output/output.nc')[source]
Run a OpenMMTools replica exchange simulation using an OpenMM coarse grained model.
- Parameters
topology – OpenMM Topology
system (System()) – OpenMM System()
positions – Positions array for the model we would like to test
total_simulation_time – Total run time for individual simulations
simulation_time_step – Simulation integration time step
temperature_list – List of temperatures for which to perform replica exchange simulations, default = None
friction – Langevin thermostat friction coefficient, default = 1 / ps
minimize (bool) – Whether minimization is done before running the simulation
output_data (netCDF4 file as generated by OpenMM) – Name of NETCDF file where we will write simulation data
exchange_frequency (int) – Number of time steps between replica exchange attempts, Default = None
output_data – file to put the output .nc
- Returns
replica_energies ( Quantity() ( np.float( [number_replicas,number_simulation_steps] ), simtk.unit ) ) - The potential energies for all replicas at all (printed) time steps
replica_positions ( Quantity() ( np.float( [number_replicas,number_simulation_steps,cgmodel.num_beads,3] ), simtk.unit ) ) - The positions for all replicas at all (printed) time steps
replica_state_indices ( np.int64( [number_replicas,number_simulation_steps] ), simtk.unit ) - The thermodynamic state assignments for all replicas at all (printed) time steps
- Example
>>> from foldamers.cg_model.cgmodel import CGModel >>> from cg_openmm.simulation.rep_exch import * >>> cgmodel = CGModel() >>> replica_energies,replica_positions,replica_state_indices = run_replica_exchange(cgmodel.topology,cgmodel.system,cgmodel.positions)