NPTSamplerMPI

class grandfep.sampler.NPTSamplerMPI(system, topology, temperature, collision_rate, timestep, log, platform=<openmm.openmm.Platform; proxy of <Swig Object of type 'OpenMM::Platform *'> >, rst_file='md.rst7', dcd_file=None, append=False, init_lambda_state=None, lambda_dict=None)

Bases: _ReplicaExchangeMixin, NPTSampler

NPT Sampler class with MPI (replica exchange) support. Only Hamiltonian is allowed to be different.

Parameters:
  • system (openmm.System) – OpenMM system object, this system should include a barostat

  • topology (app.Topology) – OpenMM topology object

  • temperature (unit.Quantity) – Reference temperature of the system, with unit

  • collision_rate (unit.Quantity) – Collision rate of the system, with unit. e.g., 1 / (1.0 * unit.picoseconds)

  • timestep (unit.Quantity) – Timestep of the simulation, with unit. e.g., 4 * unit.femtoseconds with Hydrogen Mass Repartitioning

  • log (Union[str, Path]) – Log file path for the simulation

  • platform (openmm.Platform) – OpenMM platform to use for the simulation. Default is ‘CUDA’.

  • rst_file (str) – Restart file path for the simulation. Default is “md.rst7”.

  • dcd_file (str) – DCD file path for the simulation. Default is None, which means no dcd output.

  • init_lambda_state (int) – Lambda state index for this replica, counting from 0

  • lambda_dict (dict) – A dictionary of mapping from global parameters to their values in all the sampling states.

check_temperature()

Check the reference temperature in the integrator and barostat. If they are not close, raise an Error

Returns:

temperature

Return type:

unit.Quantity

load_rst(rst_input)

Load positions/boxVector/velocities from a restart file

Parameters:

rst_input (Union[str, Path]) – Amber inpcrd file

Return type:

None

replica_exchange(calc_neighbor_only=True)

Perform one neighbor swap replica exchange. In odd RE steps, attempt exchange between 0-1, 2-3, … In even RE steps, attempt exchange between 1-2, 3-4, … If RE is accepted, update the position/boxVector/velocity

Parameters:

calc_neighbor_only (bool) – If True, only the nearest neighbor will be calculated.

Returns:

  • reduced_energy_matrix (np.array) – A numpy array with the size of (n_lambda_states, n_simulated_states)

  • re_decision (dict) – The exchange decision between all the pairs. e.g., {(0,1):(True, ratio_0_1), (2,3):(True, ratio_2_3)}. The key is the MPI rank pairs, and the value is the decision and the acceptance ratio.

replica_exchange_global_param(calc_neighbor_only=True)

Perform one neighbor swap replica exchange. In odd RE steps, attempt exchange between 0-1, 2-3, … In even RE steps, attempt exchange between 1-2, 3-4, … If RE is accepted, update the position/boxVector/velocity

Parameters:

calc_neighbor_only (bool) – If True, only the nearest neighbor will be calculated.

Returns:

  • reduced_energy_matrix (np.array) – A numpy array with the size of (n_lambda_states, n_simulated_states)

  • re_decision (dict) – The exchange decision between all the pairs. e.g., {(0,1):(True, ratio_0_1), (2,3):(True, ratio_2_3)}. The key is the MPI rank pairs, and the value is the decision and the acceptance ratio.

report_dcd(state=None)

gather the coordinates from all MPI rank, only rank 0 writes the file. Report both dcd and rst7 files.

Parameters:

state (openmm.openmm.State) – State of the simulation. If None, it will get the current state from the simulation context.

report_dcd_rank0(state=None)

gather the coordinates from all MPI rank, only rank 0 writes the file.

Parameters:

state (openmm.openmm.State) – State of the simulation. If None, it will get the current state from the simulation context.

report_rst(state=None)

gather the coordinates from all MPI rank, only rank 0 write the file.

Parameters:

state (openmm.openmm.State) – XXX

report_rst_rank0(state=None)

gather the coordinates from all MPI rank, only rank 0 write the file.

Parameters:

state (openmm.openmm.State) – XXX

set_lambda_dict(init_lambda_state, lambda_dict)
Set internal attributes :
  • init_lambda_state: (int) Lambda state index for this replica, counting from 0

  • lambda_dict: (dict) Global parameter values for all the sampling states

  • lambda_states_list: (list) Index of the Lambda state which is simulated. All the init_lambda_state in this MPI run

  • n_lambda_states: (int) Number of Lambda state to be sampled

You can have 10 states defined in the lambda_dict, and only simulate 4 replicas. In this case, each value in the lambda_dict will have a length of 10, and n_lambda_states=10, and the lambda_states_list will be a list of 4 integers, continuously increasing between 0 and 9, and init_lambda_state should be one of the values in lambda_states_list.

Parameters:
  • init_lambda_state (int) – The lambda state for this replica, counting from 0

  • lambda_dict (dict) – A dictionary of mapping from global parameters to their values in all the sampling states. For example, {"lambda_gc_vdw": [0.0, 0.5, 1.0, 1.0, 1.0], "lambda_gc_coulomb": [0.0, 0.0, 1.0, 0.5, 1.0]}

Return type:

None

set_lambda_state(i)

Set the lambda state index for this replica.

Parameters:

i (int) – The lambda state index for this replica, counting from 0

Return type:

None

set_re_step(re_step)
Parameters:

re_step (int) – Number of replica exchanges to perform

Return type:

None

set_re_step_from_log(log_input)

Read the previous log file and determine the RE step.

Parameters:

log_input (Union[str, Path]) – previous log file where the RE step counting should continue.

Return type:

None

comm

MPI communicator

dcd_reporter_dict

A dictionary of all the dcd reporter. Call the reporter inside to write the dcd trajectory file.

kBT

kB* T, with unit.

lambda_dict: dict

A dictionary of mapping from global parameters to their values in all the sampling states.

lambda_state_index: int

Lambda state index for this replica, counting from 0

lambda_states_list: list

Index of the Lambda state which is simulated. All the init_lambda_state in this MPI run can be fetched with self.lambda_states_list[rank]

logger

Logger for the Sampler

n_lambda_states: int

Number of Lambda state to be sampled. It should be equal to the length of the values in lambda_dict.

rank

Rank of this process

re_step

Number of replica exchanges performed

rst_reporter_dict

A dictionary of all the rst7 reporter. Call the reporter inside to write the rst7 restart file.

simulation: openmm.app.simulation.Simulation

Simulation ties together Topology, System, Integrator, and Context in this sampler.

size

Size of the communicator (number of processes)

system

The OpenMM System object.

temperature: openmm.unit.quantity.Quantity

reference temperature of the system, with unit

topology

The OpenMM Topology object. All the res_name, atom_index, atom_name, etc. are in this topology.