openmmtools.testsystems.LennardJonesFluid¶
- class openmmtools.testsystems.LennardJonesFluid(nparticles=1000, reduced_density=0.05, mass=Quantity(value=39.9, unit=dalton), sigma=Quantity(value=3.4, unit=angstrom), epsilon=Quantity(value=0.238, unit=kilocalorie / mole), cutoff=None, switch_width=Quantity(value=3.4, unit=angstrom), shift=False, dispersion_correction=True, lattice=False, charge=None, ewaldErrorTolerance=None, **kwargs)[source]¶
Create a periodic fluid of Lennard-Jones particles. Initial positions are assigned using a subrandom grid to minimize steric interactions.
- Parameters:
- nparticlesint, optional, default=1000
Number of Lennard-Jones particles.
- reduced_densityfloat, optional, default=0.05
Reduced density (density * sigma**3); default is appropriate for gas
- massopenmm.unit.Quantity, optional, default=39.9 * unit.amu
mass of each particle; default is appropriate for argon
- sigmaopenmm.unit.Quantity, optional, default=3.4 * unit.angstrom
Lennard-Jones sigma parameter; default is appropriate for argon
- epsilonopenmm.unit.Quantity, optional, default=0.238 * unit.kilocalories_per_mole
Lennard-Jones well depth; default is appropriate for argon
- cutoffopenmm.unit.Quantity, optional, default=None
Cutoff for nonbonded interactions. If None, defaults to 3.0 * sigma
- switch_widthopenmm.unit.Quantity with units compatible with angstroms, optional, default=3.4 * unit.angstrom
switching function is turned on at cutoff - switch_width If None, no switch will be applied (e.g. hard cutoff). Ignored if shift=True.
- shiftbool, optional, default=False
If True, will shift Lennard-Jones potential so energy will be continuous at cutoff (switch_width is ignored).
- dispersion_correctionbool, optional, default=True
if True, will use analytical dispersion correction (if not using switching function)
- latticebool, optional, default=False
If True, use fcc sphere packing to generate initial positions. The box size will be determined by nparticles and reduced_density.
- chargeopenmm.unit, optional, default=None
If not None, use alternating plus and minus charge for the particle charges. Also, if not None, use PME for electrostatics. Obviously this is no longer a traditional LJ system, but this option could be useful for testing the effect of charges in small systems.
- ewaldErrorTolerancefloat, optional, default=DEFAULT_EWALD_ERROR_TOLERANCE
The Ewald or PME tolerance. Used only if charge is not None.
Examples
Create default-size Lennard-Jones fluid.
>>> fluid = LennardJonesFluid() >>> system, positions = fluid.system, fluid.positions
Create a larger box of Lennard-Jones particles with specified reduced density.
>>> fluid = LennardJonesFluid(nparticles=1000, reduced_density=0.50) >>> system, positions = fluid.system, fluid.positions
Create Lennard-Jones fluid using switched particle interactions (switched off betwee 7 and 9 A) and more particles.
>>> fluid = LennardJonesFluid(switch_width=2.0*unit.angstroms, cutoff=9.0*unit.angstroms) >>> system, positions = fluid.system, fluid.positions
Create Lennard-Jones fluid using shifted potential.
>>> fluid = LennardJonesFluid(cutoff=9.0*unit.angstroms, shift=True) >>> system, positions = fluid.system, fluid.positions
- 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.
Methods
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__(nparticles=1000, reduced_density=0.05, mass=Quantity(value=39.9, unit=dalton), sigma=Quantity(value=3.4, unit=angstrom), epsilon=Quantity(value=0.238, unit=kilocalorie / mole), cutoff=None, switch_width=Quantity(value=3.4, unit=angstrom), shift=False, dispersion_correction=True, lattice=False, charge=None, ewaldErrorTolerance=None, **kwargs)[source]¶
Abstract base class for test system.
- Parameters:
Methods
__init__
([nparticles, reduced_density, ...])Abstract base class for test system.
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.