openmmtools.testsystems.HarmonicOscillator

class openmmtools.testsystems.HarmonicOscillator(K=Quantity(value=100.0, unit=kilocalorie / angstrom**2 * mole), mass=Quantity(value=39.948, unit=dalton), U0=Quantity(value=0.0, unit=kilojoule / mole), **kwargs)[source]

Create a 3D harmonic oscillator, with a single particle confined in an isotropic harmonic well.

Parameters:
Kopenmm.unit.Quantity, optional, default=100.0 * unit.kilocalories_per_mole/unit.angstrom**2

harmonic restraining potential

massopenmm.unit.Quantity, optional, default=39.948 * unit.amu

particle mass

U0openmm.unit.Quantity, optional, default=0.0 * unit.kilocalories_per_mole

Potential offset for harmonic oscillator

The functional form is given by
U(x) = (K/2) * ( (x-x0)^2 + y^2 + z^2 ) + U0

Notes

The natural period of a harmonic oscillator is T = 2*pi*sqrt(m/K), so you will want to use an integration timestep smaller than ~ T/10.

The standard deviation in position in each dimension is sigma = (kT / K)^(1/2)

The expectation and standard deviation of the potential energy of a 3D harmonic oscillator is (3/2)kT.

Examples

Create a 3D harmonic oscillator with default parameters:

>>> ho = HarmonicOscillator()
>>> (system, positions) = ho.system, ho.positions

Create a harmonic oscillator with specified mass and spring constant:

>>> mass = 12.0 * unit.amu
>>> K = 1.0 * unit.kilocalories_per_mole / unit.angstroms**2
>>> ho = HarmonicOscillator(K=K, mass=mass)
>>> (system, positions) = ho.system, ho.positions

Get a list of the available analytically-computed properties.

>>> print(ho.analytical_properties)
['potential_expectation', 'potential_standard_deviation']

Compute the potential expectation and standard deviation

>>> import openmm.unit as u
>>> thermodynamic_state = ThermodynamicState(temperature=298.0*u.kelvin, system=system)
>>> potential_mean = ho.get_potential_expectation(thermodynamic_state)
>>> potential_stddev = ho.get_potential_standard_deviation(thermodynamic_state)

TODO: * Add getters and setters for K, x0, U0 that access current global parameter in system * Add method to compute free energy of the harmonic oscillator(s)

Attributes:
systemopenmm.System

The openmm.System object corresponding to the test system.

positionslist

The openmm.unit.Quantity object containing the particle positions, with units compatible with openmm.unit.nanometers.

Methods

get_potential_expectation(state)

Return the expectation of the potential energy, computed analytically or numerically.

get_potential_standard_deviation(state)

Return the standard deviation of the potential energy, computed analytically or numerically.

reduced_potential_expectation(...)

Calculate the expected potential energy in state_sampled_from, divided by kB * T in state_evaluated_in.

serialize()

Return the System and positions in serialized XML form.

__init__(K=Quantity(value=100.0, unit=kilocalorie / angstrom**2 * mole), mass=Quantity(value=39.948, unit=dalton), U0=Quantity(value=0.0, unit=kilojoule / mole), **kwargs)[source]

Abstract base class for test system.

Parameters:

Methods

__init__([K, mass, U0])

Abstract base class for test system.

get_potential_expectation(state)

Return the expectation of the potential energy, computed analytically or numerically.

get_potential_standard_deviation(state)

Return the standard deviation of the potential energy, computed analytically or numerically.

reduced_potential_expectation(...)

Calculate the expected potential energy in state_sampled_from, divided by kB * T in state_evaluated_in.

serialize()

Return the System and positions in serialized XML form.

Attributes

analytical_properties

A list of available analytical properties, accessible via 'get_propertyname(thermodynamic_state)' calls.

mdtraj_topology

The mdtraj.Topology object corresponding to the test system (read-only).

name

The name of the test system.

positions

The openmm.unit.Quantity object containing the particle positions, with units compatible with openmm.unit.nanometers.

system

The openmm.System object corresponding to the test system.

topology

The openmm.app.Topology object corresponding to the test system.