openmmtools.alchemy.AlchemicalState

class openmmtools.alchemy.AlchemicalState(parameters_name_suffix=None, **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:
parameters_name_suffixstr, optional

If specified, the state will control a modified version of the global parameters with the name parameter_name + '_' + parameters_name_suffix. When this is the case, the normal parameters are not accessible.

lambda_stericsfloat, optional

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

lambda_electrostaticsfloat, optional

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

lambda_bondsfloat, optional

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

lambda_anglesfloat, optional

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

lambda_torsionsfloat, optional

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

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):
...
openmmtools.alchemy.AlchemicalStateError: The system parameter lambda_electrostatics 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.

>>> import openmm
>>> from openmm import 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

Methods

GlobalParameter(parameter_name, standard_value)

Descriptor for a global parameter.

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, *args, **kwargs)

Constructor reading the state from an alchemical system.

get_alchemical_variable(variable_name)

Return the value of the alchemical parameter.

get_function_variable(variable_name)

Return the value of the function variable.

set_alchemical_parameters(new_value)

Set all defined lambda parameters to the given value.

set_alchemical_variable(variable_name, new_value)

Set the value of the alchemical variable.

set_function_variable(variable_name, new_value)

Set the value of the function variable.

__init__(parameters_name_suffix=None, **kwargs)

Methods

__init__([parameters_name_suffix])

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, *args, **kwargs)

Constructor reading the state from an alchemical system.

get_alchemical_variable(variable_name)

Return the value of the alchemical parameter.

get_function_variable(variable_name)

Return the value of the function variable.

set_alchemical_parameters(new_value)

Set all defined lambda parameters to the given value.

set_alchemical_variable(variable_name, new_value)

Set the value of the alchemical variable.

set_function_variable(variable_name, new_value)

Set the value of the function variable.

Attributes