openmmtools.alchemy.AlchemicalState

class openmmtools.alchemy.AlchemicalState(**kwargs)[source]

Represent an alchemical state.

The alchemical parameters modify the Hamiltonian and affect the computation of the energy. Alchemical parameters that have value None are considered undefined, which means that applying this state to System and Context that have that parameter as a global variable will raise an AlchemicalStateError.

Parameters:
lambda_sterics : float, optional

Scaling factor for ligand sterics (Lennard-Jones and Halgren) interactions (default is 1.0).

lambda_electrostatics : float, optional

Scaling factor for ligand charges, intrinsic Born radii, and surface area term (default is 1.0).

lambda_bonds : float, optional

Scaling factor for alchemically-softened bonds (default is 1.0).

lambda_angles : float, optional

Scaling factor for alchemically-softened angles (default is 1.0).

lambda_torsions : float, optional

Scaling factor for alchemically-softened torsions (default is 1.0).

update_alchemical_charges : bool, optional

If False, lambda_electrostatics changes in alchemical systems that use exact treatment of PME electrostatics will be considered incompatible. This means that a new Context will be required for each lambda_electrostatics` state.

Examples

Create an alchemically modified system.

>>> from openmmtools import testsystems
>>> factory = AbsoluteAlchemicalFactory(consistent_exceptions=False)
>>> alanine_vacuum = testsystems.AlanineDipeptideVacuum().system
>>> alchemical_region = AlchemicalRegion(alchemical_atoms=range(22))
>>> alanine_alchemical_system = factory.create_alchemical_system(reference_system=alanine_vacuum,
...                                                              alchemical_regions=alchemical_region)

Create a completely undefined alchemical state.

>>> alchemical_state = AlchemicalState()
>>> print(alchemical_state.lambda_sterics)
None
>>> alchemical_state.apply_to_system(alanine_alchemical_system)
Traceback (most recent call last):
...
AlchemicalStateError: The system parameter lambda_sterics is not defined in this state.

Create an AlchemicalState that matches the parameters defined in the System.

>>> alchemical_state = AlchemicalState.from_system(alanine_alchemical_system)
>>> alchemical_state.lambda_sterics
1.0
>>> alchemical_state.lambda_electrostatics
1.0
>>> print(alchemical_state.lambda_angles)
None

AlchemicalState implement the IComposableState interface, so it can be used with CompoundThermodynamicState. All the alchemical parameters are accessible through the compound state.

>>> from simtk import openmm, unit
>>> thermodynamic_state = states.ThermodynamicState(system=alanine_alchemical_system,
...                                                 temperature=300*unit.kelvin)
>>> compound_state = states.CompoundThermodynamicState(thermodynamic_state=thermodynamic_state,
...                                                    composable_states=[alchemical_state])
>>> compound_state.lambda_sterics
1.0

You can control the parameters in the OpenMM Context in this state by setting the state attributes.

>>> compound_state.lambda_sterics = 0.5
>>> integrator = openmm.VerletIntegrator(1.0*unit.femtosecond)
>>> context = compound_state.create_context(integrator)
>>> context.getParameter('lambda_sterics')
0.5
>>> compound_state.lambda_sterics = 1.0
>>> compound_state.apply_to_context(context)
>>> context.getParameter('lambda_sterics')
1.0

You can express the alchemical parameters as a mathematical expression involving alchemical variables. Here is an example for a two-stage function.

>>> compound_state.set_alchemical_variable('lambda', 1.0)
>>> compound_state.lambda_sterics = AlchemicalFunction('step_hm(lambda - 0.5) + 2*lambda * step_hm(0.5 - lambda)')
>>> compound_state.lambda_electrostatics = AlchemicalFunction('2*(lambda - 0.5) * step(lambda - 0.5)')
>>> for l in [0.0, 0.25, 0.5, 0.75, 1.0]:
...     compound_state.set_alchemical_variable('lambda', l)
...     print(compound_state.lambda_sterics)
0.0
0.5
1.0
1.0
1.0
Attributes:
lambda_sterics
lambda_electrostatics
lambda_bonds
lambda_angles
lambda_torsions
update_alchemical_charges

Methods

apply_to_context(context) Put the Context into this AlchemicalState.
apply_to_system(system) Set the alchemical state of the system to this.
check_system_consistency(system) Check if the system is in this alchemical state.
from_system(system) Constructor reading the state from an alchemical system.
get_alchemical_variable(variable_name) Return the value of the alchemical parameter.
set_alchemical_parameters(new_value) Set all defined parameters to the given value.
set_alchemical_variable(variable_name, new_value) Set the value of the alchemical variable.
__init__(**kwargs)[source]

x.__init__(…) initializes x; see help(type(x)) for signature

Methods

__init__(**kwargs) x.__init__(…) initializes x; see help(type(x)) for signature
apply_to_context(context) Put the Context into this AlchemicalState.
apply_to_system(system) Set the alchemical state of the system to this.
check_system_consistency(system) Check if the system is in this alchemical state.
from_system(system) Constructor reading the state from an alchemical system.
get_alchemical_variable(variable_name) Return the value of the alchemical parameter.
set_alchemical_parameters(new_value) Set all defined parameters to the given value.
set_alchemical_variable(variable_name, new_value) Set the value of the alchemical variable.

Attributes

AlchemicalState.lambda_angles
AlchemicalState.lambda_bonds
AlchemicalState.lambda_electrostatics
AlchemicalState.lambda_sterics
AlchemicalState.lambda_torsions
apply_to_context(context)[source]

Put the Context into this AlchemicalState.

Parameters:
context : simtk.openmm.Context

The context to set.

Raises:
AlchemicalStateError

If the context does not have the required lambda global variables.

apply_to_system(system)[source]

Set the alchemical state of the system to this.

Parameters:
system : simtk.openmm.System

The system to modify.

Raises:
AlchemicalStateError

If the system does not have the required lambda global variables.

check_system_consistency(system)[source]

Check if the system is in this alchemical state.

It raises a AlchemicalStateError if the system is not consistent with the alchemical state.

Parameters:
system : simtk.openmm.System

The system to test.

Raises:
AlchemicalStateError

If the system is not consistent with this state.

classmethod from_system(system)[source]

Constructor reading the state from an alchemical system.

Parameters:
system : simtk.openmm.System

An alchemically modified system in a defined alchemical state.

Returns:
The AlchemicalState object representing the alchemical state of
the system.
Raises:
AlchemicalStateError

If the same parameter has different values in the system, or if the system has no lambda parameters.

get_alchemical_variable(variable_name)[source]

Return the value of the alchemical parameter.

Parameters:
variable_name : str

The name of the alchemical variable.

Returns:
variable_value : float

The value of the alchemical variable.

set_alchemical_parameters(new_value)[source]

Set all defined parameters to the given value.

The undefined parameters (i.e. those being set to None) remain undefined.

Parameters:
new_value : float

The new value for all defined parameters.

set_alchemical_variable(variable_name, new_value)[source]

Set the value of the alchemical variable.

Parameters:
variable_name : str

The name of the alchemical variable.

new_value : float

The new value for the variable.