openmmtools.states.ThermodynamicState

class openmmtools.states.ThermodynamicState(system, temperature=None, pressure=None, surface_tension=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.

NVT, NPT and NPgammaT 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.

To set a ThermodynamicState up in the NPgammaT ensemble, the system passed to the constructor has to have a MonteCarloMembraneBarostat.

To set a ThermodynamicState up with anisotropic pressure control, the system passed to the constructor has to have a MonteCarloAnisotropicBarostat. Currently the MonteCarloAnisotropicBarostat is only supported if the pressure is equal for all axes that are under pressure control.

Parameters:
systemopenmm.System

An OpenMM system in a particular thermodynamic state.

temperatureopenmm.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.

pressureopenmm.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.

surface_tensionopenmm.unit.Quantity, optional

The surface tension for the system at constant surface tension. If this is specified, the system must have a MonteCarloMembraneBarostat. If None, the surface_tension is inferred from the barostat and NPT/NVT ensemble is assumed if there is no MonteCarloMembraneBarostat.

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.

surface_tension

Surface tension

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, surface_tension=None)[source]

Methods

__init__(system[, temperature, pressure, ...])

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.

surface_tension

Surface tension

system

The system in this thermodynamic state.

temperature

Constant temperature of the thermodynamic state.

volume

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