OpenMM Relative Free Energy Protocol#

This section provides details about the OpenMM Relative Free Energy Protocol implemented in OpenFE.

Protocol Settings#

Below are the settings which can be tweaked in the protocol. The default settings (accessed using RelativeHybridTopologyProtocol.default_settings()) will automatically populate a settings which we have found to be useful for running relative binding free energies using explicit solvent. There will however be some cases (such as when doing gas phase calculations) where you will need to tweak some of the following settings.

pydantic model openfe.protocols.openmm_rfe.RelativeHybridTopologyProtocolSettings#

Create a new model by parsing and validating input data from keyword arguments.

Raises ValidationError if the input data cannot be parsed to form a valid model.

field forcefield_settings: OpenMMSystemGeneratorFFSettings [Required]#

Parameters to set up the force field with OpenMM Force Fields.

field thermo_settings: ThermoSettings [Required]#

Settings for thermodynamic parameters.

field system_settings: SystemSettings [Required]#

Simulation system settings including the long-range non-bonded method.

field solvation_settings: SolvationSettings [Required]#

Settings for solvating the system.

field alchemical_settings: AlchemicalSettings [Required]#

Alchemical protocol settings including lambda windows and soft core scaling.

field alchemical_sampler_settings: AlchemicalSamplerSettings [Required]#

Settings for sampling alchemical space, including the number of repeats.

field engine_settings: OpenMMEngineSettings [Required]#

Settings specific to the OpenMM engine such as the compute platform.

field integrator_settings: IntegratorSettings [Required]#

Settings for the integrator such as timestep and barostat settings.

field simulation_settings: SimulationSettings [Required]#

Simulation control settings, including simulation lengths and record-keeping.

pydantic model openfe.protocols.openmm_rfe.equil_rfe_settings.OpenMMSystemGeneratorFFSettings#

Parameters to set up the force field with OpenMM ForceFields

Note

Right now we just basically just grab what we need for the openmmforcefields.system_generators.SystemGenerator signature. See the OpenMMForceField SystemGenerator documentation for more details.

Create a new model by parsing and validating input data from keyword arguments.

Raises ValidationError if the input data cannot be parsed to form a valid model.

field constraints: str | None = 'hbonds'#

Constraints to be applied to system. One of ‘hbonds’, ‘allbonds’, ‘hangles’ or None, default ‘hbonds’

field rigid_water: bool = True#
field remove_com: bool = False#
field hydrogen_mass: float = 3.0#

Mass to be repartitioned to hydrogens from neighbouring heavy atoms (in amu), default 3.0

field forcefields: list[str] = ['amber/ff14SB.xml', 'amber/tip3p_standard.xml', 'amber/tip3p_HFE_multivalent.xml', 'amber/phosaa10.xml']#

List of force field paths for all components except SmallMoleculeComponent

field small_molecule_forcefield: str = 'openff-2.0.0'#

Name of the force field to be used for SmallMoleculeComponent

pydantic model openfe.protocols.openmm_rfe.equil_rfe_settings.ThermoSettings#

Settings for thermodynamic parameters.

Note

No checking is done to ensure a valid thermodynamic ensemble is possible.

Create a new model by parsing and validating input data from keyword arguments.

Raises ValidationError if the input data cannot be parsed to form a valid model.

field temperature: FloatQuantity = None#

Simulation temperature, default units kelvin

field pressure: FloatQuantity = None#

Simulation pressure, default units standard atmosphere (atm)

field ph: PositiveFloat | None = None#

Simulation pH

Constraints:
  • exclusiveMinimum = 0

field redox_potential: float | None = None#

Simulation redox potential

pydantic model openfe.protocols.openmm_rfe.equil_rfe_settings.AlchemicalSamplerSettings#

Settings for the Equilibrium Alchemical sampler, currently supporting either MultistateSampler, SAMSSampler or ReplicaExchangeSampler.

Create a new model by parsing and validating input data from keyword arguments.

Raises ValidationError if the input data cannot be parsed to form a valid model.

field sampler_method = 'repex'#

Alchemical sampling method, must be one of; repex (Hamiltonian Replica Exchange), sams (Self-Adjusted Mixture Sampling), or independent (independently sampled lambda windows). Default repex.

field online_analysis_interval: int | None = 250#

MCMC steps (i.e. IntegratorSettings.n_steps) interval at which to perform an analysis of the free energies.

At each interval, real time analysis data (e.g. current free energy estimate and timing data) will be written to a yaml file named <SimulationSettings.output_name>_real_time_analysis.yaml. The current error in the estimate will also be assed and if it drops below AlchemicalSamplerSettings.online_analysis_target_error the simulation will be terminated.

If None, no real time analysis will be performed and the yaml file will not be written.

Must be a multiple of SimulationSettings.checkpoint_interval

Default 250.

field online_analysis_target_error = <Quantity(0.0, 'boltzmann_constant * kelvin')>#

Target error for the online analysis measured in kT. Once the free energy is at or below this value, the simulation will be considered complete. A suggested value of 0.2 * unit.boltzmann_constant * unit.kelvin has shown to be effective in both hydration and binding free energy benchmarks. Default 0.0 * unit.boltzmann_constant * unit.kelvin, i.e. no early termination will occur.

field online_analysis_minimum_iterations = 500#

Number of iterations which must pass before online analysis is carried out. Default 500.

field n_repeats: int = 3#

Number of independent repeats to run. Default 3

field flatness_criteria = 'logZ-flatness'#

SAMS only. Method for assessing when to switch to asymptomatically optimal scheme. One of [‘logZ-flatness’, ‘minimum-visits’, ‘histogram-flatness’]. Default ‘logZ-flatness’.

field gamma0 = 1.0#

SAMS only. Initial weight adaptation rate. Default 1.0.

field n_replicas = 11#

Number of replicas to use. Default 11.

pydantic model openfe.protocols.openmm_rfe.equil_rfe_settings.AlchemicalSettings#

Settings for the alchemical protocol

This describes the lambda schedule and the creation of the hybrid system.

Create a new model by parsing and validating input data from keyword arguments.

Raises ValidationError if the input data cannot be parsed to form a valid model.

field lambda_functions = 'default'#

Key of which switching functions to use for alchemical mutation. Default ‘default’.

field lambda_windows = 11#

Number of lambda windows to calculate. Default 11.

field unsampled_endstates = False#

Whether to have extra unsampled endstate windows for long range correction. Default False.

field use_dispersion_correction = False#

Whether to use dispersion correction in the hybrid topology state. Default False.

field softcore_LJ_v2 = True#

Whether to use the LJ softcore function as defined by Gapsys et al. JCTC 2012 Default True.

field softcore_electrostatics = True#

Whether to use softcore electrostatics. Default True.

field softcore_alpha = 0.85#

Softcore alpha parameter. Default 0.85

field softcore_electrostatics_alpha = 0.3#

Softcore alpha parameter for electrostatics. Default 0.3

field softcore_sigma_Q = 1.0#

Softcore sigma parameter for softcore electrostatics. Default 1.0.

field interpolate_old_and_new_14s = False#

Whether to turn off interactions for new exceptions (not just 1,4s) at lambda 0 and old exceptions at lambda 1. If False they are present in the nonbonded force. Default False.

field flatten_torsions = False#

Whether to scale torsion terms involving unique atoms, such that at the endstate the torsion term is turned off/on depending on the state the unique atoms belong to. Default False.

pydantic model openfe.protocols.openmm_rfe.equil_rfe_settings.OpenMMEngineSettings#

OpenMM MD engine settings

Create a new model by parsing and validating input data from keyword arguments.

Raises ValidationError if the input data cannot be parsed to form a valid model.

field compute_platform: str | None = None#

OpenMM compute platform to perform MD integration with. If None, will choose fastest available platform. Default None.

pydantic model openfe.protocols.openmm_rfe.equil_rfe_settings.IntegratorSettings#

Settings for the LangevinSplittingDynamicsMove integrator

Create a new model by parsing and validating input data from keyword arguments.

Raises ValidationError if the input data cannot be parsed to form a valid model.

field timestep = <Quantity(4, 'femtosecond')>#

Size of the simulation timestep. Default 4 * unit.femtosecond.

field collision_rate = <Quantity(1.0, '1 / picosecond')>#

Collision frequency. Default 1.0 / unit.pisecond.

field n_steps = <Quantity(250, 'timestep')>#

Number of integration timesteps between each time the MCMC move is applied. Default 250 * unit.timestep.

field reassign_velocities = False#

If True, velocities are reassigned from the Maxwell-Boltzmann distribution at the beginning of move. Default False.

field n_restart_attempts = 20#

Number of attempts to restart from Context if there are NaNs in the energies after integration. Default 20.

field constraint_tolerance = 1e-06#

Tolerance for the constraint solver. Default 1e-6.

field barostat_frequency = <Quantity(25, 'timestep')>#

Frequency at which volume scaling changes should be attempted. Default 25 * unit.timestep.

pydantic model openfe.protocols.openmm_rfe.equil_rfe_settings.SimulationSettings#

Settings for simulation control, including lengths, writing to disk, etc…

Create a new model by parsing and validating input data from keyword arguments.

Raises ValidationError if the input data cannot be parsed to form a valid model.

field minimization_steps = 5000#

Number of minimization steps to perform. Default 5000.

field equilibration_length: Quantity [Required]#

Length of the equilibration phase in units of time. The total number of steps from this equilibration length (i.e. equilibration_length / IntegratorSettings.timestep) must be a multiple of the value defined for IntegratorSettings.n_steps.

field production_length: Quantity [Required]#

Length of the production phase in units of time. The total number of steps from this production length (i.e. production_length / IntegratorSettings.timestep) must be a multiple of the value defined for IntegratorSettings.nsteps.

field output_filename = 'simulation.nc'#

Path to the trajectory storage file. Default ‘simulation.nc’.

field output_structure = 'hybrid_system.pdb'#

Path of the output hybrid topology structure file. This is used to visualise and further manipulate the system. Default ‘hybrid_system.pdb’.

field output_indices = 'not water'#

Selection string for which part of the system to write coordinates for. Default ‘not water’.

field checkpoint_interval = <Quantity(250, 'timestep')>#

Frequency to write the checkpoint file. Default 250 * unit.timestep.

field checkpoint_storage = 'checkpoint.nc'#

Separate filename for the checkpoint file. Note, this should not be a full path, just a filename. Default ‘checkpoint.nc’.

field forcefield_cache: str | None = 'db.json'#

Filename for caching small molecule residue templates so they can be later reused.

pydantic model openfe.protocols.openmm_rfe.equil_rfe_settings.SolvationSettings#

Settings for solvating the system

Note

No solvation will happen if a SolventComponent is not passed.

Create a new model by parsing and validating input data from keyword arguments.

Raises ValidationError if the input data cannot be parsed to form a valid model.

field solvent_model = 'tip3p'#

Force field water model to use. Allowed values are; tip3p, spce, tip4pew, and tip5p.

field solvent_padding = <Quantity(1.2, 'nanometer')>#

Minimum distance from any solute atoms to the solvent box edge.

pydantic model openfe.protocols.openmm_rfe.equil_rfe_settings.SystemSettings#

Settings describing the simulation system settings.

Create a new model by parsing and validating input data from keyword arguments.

Raises ValidationError if the input data cannot be parsed to form a valid model.

field nonbonded_method = 'PME'#

Method for treating nonbonded interactions, currently only PME and NoCutoff are allowed. Default PME.

field nonbonded_cutoff = <Quantity(1.0, 'nanometer')>#

Cutoff value for short range nonbonded interactions. Default 1.0 * unit.nanometer.

Protocol API specification#

class openfe.protocols.openmm_rfe.RelativeHybridTopologyProtocol(settings: Settings)#

Create a new Protocol instance.

Parameters:

settings (Settings) – The full settings for this Protocol instance.

classmethod default_settings() Settings#

Get the default settings for this Protocol.

These can be modified and passed in as the settings for a new Protocol instance.

create(*, stateA: ChemicalSystem, stateB: ChemicalSystem, mapping: dict[str, ComponentMapping] | None, extends: ProtocolDAGResult | None = None, name: str | None = None, transformation_key: GufeKey | None = None) ProtocolDAG#

Prepare a ProtocolDAG with all information required for execution.

A ProtocolDAG is composed of ProtocolUnit`s, with dependencies established between them. These form a directed, acyclic graph, and each `ProtocolUnit can be executed once its dependencies have completed.

A ProtocolDAG can be passed to a Scheduler for execution on its resources. A ProtocolDAGResult can be retrieved from the Scheduler upon completion of all ProtocolUnit`s in the `ProtocolDAG.

Parameters:
  • stateA (ChemicalSystem) – The starting ChemicalSystem for the transformation.

  • stateB (ChemicalSystem) – The ending ChemicalSystem for the transformation.

  • mapping (Optional[dict[str, ComponentMapping]]) –

    Mappings of e.g. atoms between a labelled component in the

    stateA and stateB ChemicalSystem .

  • extends (Optional[ProtocolDAGResult]) – If provided, then the ProtocolDAG produced will start from the end state of the given ProtocolDAGResult. This allows for extension from a previously-run ProtocolDAG.

  • name (Optional[str]) – A user supplied identifier for the resulting DAG

  • transformation_key (Optional[GufeKey]) – Key of the Transformation that this Protocol corresponds to, if applicable. This will be used to label the resulting ProtocolDAG, and can be used for identifying its source. This label will be passed on to the ProtocolDAGResult resulting from execution of this ProtocolDAG.

Returns:

A directed, acyclic graph that can be executed by a Scheduler.

Return type:

ProtocolDAG

gather(protocol_dag_results: Iterable[ProtocolDAGResult]) ProtocolResult#

Gather multiple ProtocolDAGResults into a single ProtocolResult.

Parameters:

protocol_dag_results (Iterable[ProtocolDAGResult]) – The ProtocolDAGResult objects to assemble aggregate quantities from.

Returns:

Aggregated results from many ProtocolDAGResult`s from a given `Protocol.

Return type:

ProtocolResult

class openfe.protocols.openmm_rfe.RelativeHybridTopologyProtocolResult(**data)#

Dict-like container for the output of a RelativeHybridTopologyProtocol

get_estimate()#

Average free energy difference of this transformation

Returns:

dG – The free energy difference between the first and last states. This is a Quantity defined with units.

Return type:

unit.Quantity

get_uncertainty()#

The uncertainty/error in the dG value: The std of the estimates of each independent repeat