openmmtools.testsystems.HarmonicOscillatorArray

class openmmtools.testsystems.HarmonicOscillatorArray(K=Quantity(value=90.0, unit=kilocalorie / angstrom**2 * mole), d=Quantity(value=1.0, unit=nanometer), mass=Quantity(value=39.948, unit=dalton), N=5, **kwargs)[source]

Create a 1D array of noninteracting particles in 3D harmonic oscillator wells.

Parameters:
Kopenmm.unit.Quantity, optional, default=90.0 * unit.kilocalories_per_mole/unit.angstroms**2

harmonic restraining potential

dopenmm.unit.Quantity, optional, default=1.0 * unit.nanometer

distance between harmonic oscillators. Default is Amber GAFF c-c bond.

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

particle mass

Nint, optional, default=5

Number of harmonic oscillators

Notes

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

Examples

Create a constraint-coupled 3D harmonic oscillator with default parameters.

>>> ho_array = HarmonicOscillatorArray()
>>> mass = 12.0 * unit.amu
>>> d = 5.0 * unit.angstroms
>>> K = 1.0 * unit.kilocalories_per_mole / unit.angstroms**2
>>> ccho = HarmonicOscillatorArray(K=K, d=d, mass=mass)
>>> system, positions = ccho.system, ccho.positions
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=90.0, unit=kilocalorie / angstrom**2 * mole), d=Quantity(value=1.0, unit=nanometer), mass=Quantity(value=39.948, unit=dalton), N=5, **kwargs)[source]

Abstract base class for test system.

Parameters:

Methods

__init__([K, d, mass, N])

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.