openmmtools.states.ThermodynamicState

class openmmtools.states.ThermodynamicState(system, temperature=None, pressure=None)[source]

Thermodynamic state of a system.

Represent the portion of the state of a Context that does not change with integration. Its main objectives are to wrap an OpenMM system object to easily maintain a consistent thermodynamic state. It can be used to create new OpenMM Contexts, or to convert an existing Context to this particular thermodynamic state.

Only NVT and NPT ensembles are supported. The temperature must be specified in the constructor, either implicitly via a thermostat force in the system, or explicitly through the temperature parameter, which overrides an eventual thermostat indication.

Parameters:
system : simtk.openmm.System

An OpenMM system in a particular thermodynamic state.

temperature : simtk.unit.Quantity, optional

The temperature for the system at constant temperature. If a MonteCarloBarostat is associated to the system, its temperature will be set to this. If None, the temperature is inferred from the system thermostat.

pressure : simtk.unit.Quantity, optional

The pressure for the system at constant pressure. If this is specified, a MonteCarloBarostat is added to the system, or just set to this pressure in case it already exists. If None, the pressure is inferred from the system barostat, and NVT ensemble is assumed if there is no barostat.

Notes

This state object cannot describe states obeying non-Boltzamnn statistics, such as Tsallis statistics.

Examples

Specify an NVT state for a water box at 298 K.

>>> from openmmtools import testsystems
>>> temperature = 298.0*unit.kelvin
>>> waterbox = testsystems.WaterBox(box_edge=10*unit.angstroms,
...                                 cutoff=4*unit.angstroms).system
>>> state = ThermodynamicState(system=waterbox, temperature=temperature)

In an NVT ensemble volume is constant and pressure is None.

>>> state.volume
Quantity(value=1.0, unit=nanometer**3)
>>> state.pressure is None
True

Convert this to an NPT state at 298 K and 1 atm pressure. This operation automatically adds a MonteCarloBarostat to the system.

>>> pressure = 1.0*unit.atmosphere
>>> state.pressure = pressure
>>> state.pressure
Quantity(value=1.0, unit=atmosphere)
>>> state.volume is None
True

You cannot set a non-periodic system at constant pressure

>>> nonperiodic_system = testsystems.TolueneVacuum().system
>>> state = ThermodynamicState(nonperiodic_system, temperature=300*unit.kelvin,
...                            pressure=1.0*unit.atmosphere)
Traceback (most recent call last):
...
openmmtools.states.ThermodynamicsError: Non-periodic systems cannot have a barostat.

When temperature and/or pressure are not specified (i.e. they are None) ThermodynamicState tries to infer them from a thermostat or a barostat.

>>> state = ThermodynamicState(system=waterbox)
Traceback (most recent call last):
...
openmmtools.states.ThermodynamicsError: System does not have a thermostat specifying the temperature.
>>> thermostat = openmm.AndersenThermostat(200.0*unit.kelvin, 1.0/unit.picosecond)
>>> force_id = waterbox.addForce(thermostat)
>>> state = ThermodynamicState(system=waterbox)
>>> state.pressure is None
True
>>> state.temperature
Quantity(value=200.0, unit=kelvin)
>>> barostat = openmm.MonteCarloBarostat(1.0*unit.atmosphere, 200.0*unit.kelvin)
>>> force_id = waterbox.addForce(barostat)
>>> state = ThermodynamicState(system=waterbox)
>>> state.pressure
Quantity(value=1.01325, unit=bar)
>>> state.temperature
Quantity(value=200.0, unit=kelvin)
Attributes:
system

The system in this thermodynamic state.

temperature

Constant temperature of the thermodynamic state.

pressure

Constant pressure of the thermodynamic state.

volume

Constant volume of the thermodynamic state (read-only).

n_particles

Number of particles (read-only).

Methods

apply_to_context(context) Apply this ThermodynamicState to the context.
create_context(integrator[, platform]) Create a context in this ThermodynamicState.
get_system([remove_thermostat, remove_barostat]) Manipulate and return the system.
get_volume([ignore_ensemble]) Volume of the periodic box (read-only).
is_context_compatible(context) Check compatibility of the given context.
is_state_compatible(thermodynamic_state) Check compatibility between ThermodynamicStates.
reduced_potential(context_state) Reduced potential in this thermodynamic state.
reduced_potential_at_states(context, …) Efficiently compute the reduced potential for a list of compatible states.
set_system(system[, fix_state]) Manipulate and set the system.
__init__(system, temperature=None, pressure=None)[source]

Initialize self. See help(type(self)) for accurate signature.

Methods

__init__(system[, temperature, pressure]) Initialize self.
apply_to_context(context) Apply this ThermodynamicState to the context.
create_context(integrator[, platform]) Create a context in this ThermodynamicState.
get_system([remove_thermostat, remove_barostat]) Manipulate and return the system.
get_volume([ignore_ensemble]) Volume of the periodic box (read-only).
is_context_compatible(context) Check compatibility of the given context.
is_state_compatible(thermodynamic_state) Check compatibility between ThermodynamicStates.
reduced_potential(context_state) Reduced potential in this thermodynamic state.
reduced_potential_at_states(context, …) Efficiently compute the reduced potential for a list of compatible states.
set_system(system[, fix_state]) Manipulate and set the system.

Attributes

barostat The barostat associated to the system.
beta Thermodynamic beta in units of mole/energy.
default_box_vectors The default box vectors of the System (read-only).
is_periodic True if the system is in a periodic box (read-only).
kT Thermal energy per mole.
n_particles Number of particles (read-only).
pressure Constant pressure of the thermodynamic state.
system The system in this thermodynamic state.
temperature Constant temperature of the thermodynamic state.
volume Constant volume of the thermodynamic state (read-only).
apply_to_context(context)[source]

Apply this ThermodynamicState to the context.

The method apply_to_context does not check for the compatibility of the context. The user is responsible for this. Depending on the system size, is_context_compatible can be an expensive operation, so is_state_compatible should be preferred when possible.

Parameters:
context : simtk.openmm.Context

The OpenMM Context to be set to this ThermodynamicState.

Raises:
ThermodynamicsError

If the context is in a different thermodynamic ensemble w.r.t. this state. This is just a quick check which does not substitute is_state_compatible or is_context_compatible.

Examples

The method doesn’t verify compatibility with the context, it is the user’s responsibility to do so, possibly with is_state_compatible rather than is_context_compatible which is slower.

>>> from simtk import openmm, unit
>>> from openmmtools import testsystems
>>> toluene = testsystems.TolueneVacuum()
>>> state1 = ThermodynamicState(toluene.system, 273.0*unit.kelvin)
>>> state2 = ThermodynamicState(toluene.system, 310.0*unit.kelvin)
>>> integrator = openmm.VerletIntegrator(1.0*unit.femtosecond)
>>> context = state1.create_context(integrator)
>>> if state2.is_state_compatible(state1):
...     state2.apply_to_context(context)
>>> context.getParameter(openmm.AndersenThermostat.Temperature())
310.0
barostat

The barostat associated to the system.

Note that this is only a copy of the barostat, and you will need to set back the ThermodynamicState.barostat property for the changes to take place internally. If the pressure is allowed to fluctuate, this is None. Normally, you should only need to access the pressure and temperature properties, but this allows you to modify other parameters of the MonteCarloBarostat (e.g. frequency) after initialization. Setting this to None will place the system in an NVT ensemble.

beta

Thermodynamic beta in units of mole/energy.

create_context(integrator, platform=None)[source]

Create a context in this ThermodynamicState.

The context contains a copy of the system. If the integrator is coupled to a heat bath (e.g. LangevinIntegrator), the system in the context will not have a thermostat, and vice versa if the integrator is not thermostated the system in the context will have a thermostat.

An integrator is considered thermostated if it exposes a method getTemperature(). A CompoundIntegrator is considered coupled to a heat bath if at least one of its integrators is. An exception is raised if the integrator is thermostated at a temperature different from the thermodynamic state’s.

Parameters:
integrator : simtk.openmm.Integrator

The integrator to use for Context creation. The eventual heat bath temperature must be consistent with the thermodynamic state.

platform : simtk.openmm.Platform, optional

Platform to use. If None, OpenMM tries to select the fastest available platform. Default is None.

Returns:
context : simtk.openmm.Context

The created OpenMM Context object.

Raises:
ThermodynamicsError

If the integrator has a temperature different from this ThermodynamicState.

Examples

When passing an integrator that does not expose getter and setter for the temperature, the context will be created with a thermostat.

>>> from simtk import openmm, unit
>>> from openmmtools import testsystems
>>> toluene = testsystems.TolueneVacuum()
>>> state = ThermodynamicState(toluene.system, 300*unit.kelvin)
>>> integrator = openmm.VerletIntegrator(1.0*unit.femtosecond)
>>> context = state.create_context(integrator)
>>> system = context.getSystem()
>>> [force.__class__.__name__ for force in system.getForces()
...  if 'Thermostat' in force.__class__.__name__]
['AndersenThermostat']

The thermostat is removed if we choose an integrator coupled to a heat bath.

>>> del context  # Delete previous context to free memory.
>>> integrator = openmm.LangevinIntegrator(300*unit.kelvin, 5.0/unit.picosecond,
...                                        2.0*unit.femtosecond)
>>> context = state.create_context(integrator)
>>> system = context.getSystem()
>>> [force.__class__.__name__ for force in system.getForces()
...  if 'Thermostat' in force.__class__.__name__]
[]
default_box_vectors

The default box vectors of the System (read-only).

get_system(remove_thermostat=False, remove_barostat=False)[source]

Manipulate and return the system.

With default arguments, this is equivalent as the system property. By setting the arguments it is possible to obtain a modified copy of the system without the thermostat or the barostat.

Parameters:
remove_thermostat : bool

If True, the system thermostat is removed.

remove_barostat : bool

If True, the system barostat is removed.

Returns:
system : simtk.openmm.System

The system of this ThermodynamicState.

Examples

The constructor adds a thermostat and a barostat to configure the system in an NPT ensemble.

>>> from openmmtools import testsystems
>>> alanine = testsystems.AlanineDipeptideExplicit()
>>> state = ThermodynamicState(alanine.system, temperature=300*unit.kelvin,
...                            pressure=1.0*unit.atmosphere)

The system property returns a copy of the system with the added thermostat and barostat.

>>> system = state.system
>>> [force.__class__.__name__ for force in system.getForces()
...  if 'Thermostat' in force.__class__.__name__]
['AndersenThermostat']

We can remove them while getting the arguments with

>>> system = state.get_system(remove_thermostat=True, remove_barostat=True)
>>> [force.__class__.__name__ for force in system.getForces()
...  if 'Thermostat' in force.__class__.__name__]
[]
get_volume(ignore_ensemble=False)[source]

Volume of the periodic box (read-only).

Parameters:
ignore_ensemble : bool, optional

If True, the volume of the periodic box vectors is returned even if the volume fluctuates.

Returns:
volume : simtk.unit.Quantity

The volume of the periodic box (units of length^3) or None if the system is not periodic or allowed to fluctuate.

is_context_compatible(context)[source]

Check compatibility of the given context.

This is equivalent to is_state_compatible but slower, and it should be used only when the state the created the context is unknown. The context is compatible if it can be set to this ThermodynamicState through apply_to_context().

Parameters:
context : simtk.openmm.Context

The OpenMM context to test.

Returns:
is_compatible : bool

True if this ThermodynamicState can be applied to context.

is_periodic

True if the system is in a periodic box (read-only).

is_state_compatible(thermodynamic_state)[source]

Check compatibility between ThermodynamicStates.

The state is compatible if Contexts created by thermodynamic_state can be set to this ThermodynamicState through apply_to_context. The property is symmetric and transitive.

This is faster than checking compatibility of a Context object through is_context_compatible, and it should be preferred when possible.

Parameters:
thermodynamic_state : ThermodynamicState

The thermodynamic state to test.

Returns:
is_compatible : bool

True if the context created by thermodynamic_state can be converted to this state through apply_to_context().

Examples

States in the same ensemble (NVT or NPT) are compatible.

>>> from simtk import unit
>>> from openmmtools import testsystems
>>> alanine = testsystems.AlanineDipeptideExplicit()
>>> state1 = ThermodynamicState(alanine.system, 273*unit.kelvin)
>>> state2 = ThermodynamicState(alanine.system, 310*unit.kelvin)
>>> state1.is_state_compatible(state2)
True

States in different ensembles are not compatible.

>>> state1.pressure = 1.0*unit.atmosphere
>>> state1.is_state_compatible(state2)
False

States that store different systems (that differ by more than barostat and thermostat pressure and temperature) are also not compatible.

>>> alanine_implicit = testsystems.AlanineDipeptideImplicit().system
>>> state_implicit = ThermodynamicState(alanine_implicit, 310*unit.kelvin)
>>> state2.is_state_compatible(state_implicit)
False
kT

Thermal energy per mole.

n_particles

Number of particles (read-only).

pressure

Constant pressure of the thermodynamic state.

If the pressure is allowed to fluctuate, this is None. Setting this will automatically add/configure a barostat to the system. If it is set to None, the barostat will be removed.

reduced_potential(context_state)[source]

Reduced potential in this thermodynamic state.

Parameters:
context_state : SamplerState or simtk.openmm.Context

Carry the configurational properties of the system.

Returns:
u : float

The unit-less reduced potential, which can be considered to have units of kT.

Notes

The reduced potential is defined as in Ref. [1]

u = eta [U(x) + p V(x) + mu N(x)]

where the thermodynamic parameters are

eta = 1/(kB T) is the inverse temperature p is the pressure mu is the chemical potential

and the configurational properties are

x the atomic positions U(x) is the potential energy V(x) is the instantaneous box volume N(x) the numbers of various particle species (e.g. protons of

titratible groups)

References

[1] Shirts MR and Chodera JD. Statistically optimal analysis of equilibrium states. J Chem Phys 129:124105, 2008.

Examples

Compute the reduced potential of a water box at 298 K and 1 atm.

>>> from openmmtools import testsystems
>>> waterbox = testsystems.WaterBox(box_edge=20.0*unit.angstroms)
>>> system, positions = waterbox.system, waterbox.positions
>>> state = ThermodynamicState(system=waterbox.system,
...                            temperature=298.0*unit.kelvin,
...                            pressure=1.0*unit.atmosphere)
>>> integrator = openmm.VerletIntegrator(1.0*unit.femtosecond)
>>> context = state.create_context(integrator)
>>> context.setPositions(waterbox.positions)
>>> sampler_state = SamplerState.from_context(context)
>>> u = state.reduced_potential(sampler_state)

If the sampler state is incompatible, an error is raised

>>> incompatible_sampler_state = sampler_state[:-1]
>>> state.reduced_potential(incompatible_sampler_state)
Traceback (most recent call last):
...
openmmtools.states.ThermodynamicsError: The sampler state has a different number of particles.

In case a cached SamplerState containing the potential energy and the volume of the context is not available, the method accepts a Context object and compute them with Context.getState().

>>> u = state.reduced_potential(context)
classmethod reduced_potential_at_states(context, thermodynamic_states)[source]

Efficiently compute the reduced potential for a list of compatible states.

The user is responsible to ensure that the given context is compatible with the thermodynamic states.

Parameters:
context : openmm.Context

The OpenMM Context object with box vectors and positions set.

thermodynamic_states : list of ThermodynamicState

The list of thermodynamic states at which to compute the reduced potential.

Returns:
reduced_potentials : list of float

The unit-less reduced potentials, which can be considered to have units of kT.

Raises:
ValueError

If the thermodynamic states are not compatible to each other.

set_system(system, fix_state=False)[source]

Manipulate and set the system.

With default arguments, this is equivalent to using the system property, which raises an exception if the thermostat and the barostat are not configured according to the thermodynamic state. With this method it is possible to adjust temperature and pressure of the system to make the assignment possible, without manually configuring thermostat and barostat.

Parameters:
system : simtk.openmm.System

The system to set.

fix_state : bool, optional

If True, a thermostat is added to the system (if not already present) and set to the correct temperature. If this state is in NPT ensemble, a barostat is added or configured if it exist already. If False, this simply check that thermostat and barostat are correctly configured without modifying them. Default is False.

Raises:
ThermodynamicsError

If the system after the requested manipulation is still in an incompatible state.

Examples

The constructor adds a thermostat and a barostat to configure the system in an NPT ensemble.

>>> from openmmtools import testsystems
>>> alanine = testsystems.AlanineDipeptideExplicit()
>>> state = ThermodynamicState(alanine.system, temperature=300*unit.kelvin,
...                            pressure=1.0*unit.atmosphere)

If we try to set a system not in NPT ensemble, an error occur.

>>> state.system = alanine.system
Traceback (most recent call last):
...
openmmtools.states.ThermodynamicsError: System does not have a thermostat specifying the temperature.

We can fix both thermostat and barostat while setting the system. >>> state.set_system(alanine.system, fix_state=True)

system

The system in this thermodynamic state.

The returned system is a copy and can be modified without altering the internal state of ThermodynamicState. In order to ensure a consistent thermodynamic state, the system has a Thermostat force. You can use get_system() to obtain a copy of the system without the thermostat. The method create_context() then takes care of removing the thermostat when an integrator with a coupled heat bath is used (e.g. LangevinIntegrator).

It can be set only to a system which is consistent with the current thermodynamic state. Use set_system() if you want to correct the thermodynamic state of the system automatically before assignment.

temperature

Constant temperature of the thermodynamic state.

volume

Constant volume of the thermodynamic state (read-only).

If the volume is allowed to fluctuate, or if the system is not in a periodic box this is None.