openmmtools.testsystems.Diatom¶
- class openmmtools.testsystems.Diatom(K=290.1 kcal/(A**2 mol), r0=1.55 A, m1=39.948 Da, m2=39.948 Da, constraint=False, use_central_potential=False, **kwargs)[source]¶
Create a free diatomic molecule with a single harmonic bond between the two atoms.
- Parameters:
- Kopenmm.unit.Quantity, optional, default=290.1 * unit.kilocalories_per_mole / unit.angstrom**2
harmonic bond potential. default is GAFF c-c bond
- r0openmm.unit.Quantity, optional, default=1.550 * unit.amu
bond length. Default is Amber GAFF c-c bond.
- constraintbool, default=False
if True, the bond length will be constrained
- m1openmm.unit.Quantity, optional, default=12.01 * unit.amu
particle1 mass
- m2openmm.unit.Quantity, optional, default=12.01 * unit.amu
particle2 mass
- use_central_potentialbool, optional, default=False
if True, a soft central potential will also be added to keep the system from drifting away
- Attributes:
analytical_propertiesA list of available analytical properties, accessible via ‘get_propertyname(thermodynamic_state)’ calls.
mdtraj_topologyThe mdtraj.Topology object corresponding to the test system (read-only).
nameThe name of the test system.
positionsThe openmm.unit.Quantity object containing the particle positions, with units compatible with openmm.unit.nanometers.
systemThe openmm.System object corresponding to the test system.
topologyThe openmm.app.Topology object corresponding to the test system.
Methods
get_potential_expectation(state)Return the expectation 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.
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 Diatom:
>>> diatom = Diatom() >>> system, positions = diatom.system, diatom.positions
Create a Diatom with constraint in a central potential >>> diatom = Diatom(constraint=True, use_central_potential=True) >>> system, positions = diatom.system, diatom.positions
- __init__(K=290.1 kcal/(A**2 mol), r0=1.55 A, m1=39.948 Da, m2=39.948 Da, constraint=False, use_central_potential=False, **kwargs)[source]¶
Abstract base class for test system.
- Parameters:
Methods
__init__([K, r0, m1, m2, constraint, ...])Abstract base class for test system.
get_potential_expectation(state)Return the expectation 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_propertiesA list of available analytical properties, accessible via 'get_propertyname(thermodynamic_state)' calls.
mdtraj_topologyThe mdtraj.Topology object corresponding to the test system (read-only).
nameThe name of the test system.
positionsThe openmm.unit.Quantity object containing the particle positions, with units compatible with openmm.unit.nanometers.
systemThe openmm.System object corresponding to the test system.
topologyThe openmm.app.Topology object corresponding to the test system.