TinkerFiles

class openmm.app.tinkerfiles.TinkerFiles(xyzFile: str, parameterFiles: str | List[str], periodicBoxVectors: Tuple[Vec3, Vec3, Vec3] | None = None, unitCellDimensions: Vec3 | None = None)

TinkerFiles parses Tinker files (.xyz, .prm, .key), constructs a Topology, and (optionally) an OpenMM System from it. This class only supports the AMOEBA force field. It cannot create a System from Tinker files that use other force fields.

__init__(xyzFile: str, parameterFiles: str | List[str], periodicBoxVectors: Tuple[Vec3, Vec3, Vec3] | None = None, unitCellDimensions: Vec3 | None = None)

Load a set of Tinker files, including one .xyz files and one or more .key or .prm files.

Parameters:
  • xyzFile (str) – Path to the Tinker .xyz file.

  • parameterFiles (str or list[str]) – Paths to the Tinker .key and .prm files.

  • periodicBoxVectors (Optional[Tuple[Vec3, Vec3, Vec3]]) – The periodic box vectors.

  • unitCellDimensions (Optional[Vec3]) – The unit cell dimensions.

Methods

__init__(xyzFile, parameterFiles[, ...])

Load a set of Tinker files, including one .xyz files and one or more .key or .prm files.

createSystem([nonbondedMethod, ...])

Create an OpenMM System from the parsed Tinker files.

getBoxVectors([asNumpy])

Get the periodic box vectors.

getPositions([asNumpy])

Get the atomic positions.

Attributes

GK_PARAMS

RECOGNIZED_FORCES

RECOGNIZED_SCALARS

WCA_PARAMS

createSystem(nonbondedMethod=NoCutoff, nonbondedCutoff=1.0 nm, vdwCutoff=None, removeCMMotion: bool = True, hydrogenMass=None, polarization: str = 'mutual', mutualInducedTargetEpsilon: float = 1e-05, implicitSolvent: bool = False, ewaldErrorTolerance=0.0005, *args, **kwargs) Any

Create an OpenMM System from the parsed Tinker files.

Parameters:
  • nonbondedMethod (object=PME) – The method to use for nonbonded interactions. Allowed values are NoCutoff, and PME.

  • nonbondedCutoff (distance=1.0*nanometers) – The cutoff distance to use for nonbonded interactions.

  • vdwCutoff (distance=None) – An optional alternate cutoff to use for vdw interactions. If this is omitted, nonbondedCutoff is used for vdw interactions as well as for other nonbonded interactions.

  • removeCMMotion (bool, optional, default=True) – If True, center of mass motion will be removed.

  • hydrogenMass (mass=None) – The mass to use for hydrogen atoms bound to heavy atoms. Any mass added to a hydrogen is subtracted from the heavy atom to keep their total mass the same.

  • polarization (str, optional, default="mutual") – The method to use for calculating induced dipoles. Allowed values are “mutual”, “direct”, or “extrapolated”.

  • mutualInducedTargetEpsilon (float, optional, default=0.00001) – The target epsilon for mutual induced dipoles. Only used if polarization=”mutual”.

  • implicitSolvent (bool, optional, default=False) – If True, solvent will be modeled implicitly.

  • ewaldErrorTolerance (float=0.0005) – The error tolerance to use if nonbondedMethod is Ewald, PME, or LJPME.

Returns:

The created OpenMM System.

Return type:

openmm.System

getPositions(asNumpy: bool = False) List[Vec3] | ndarray

Get the atomic positions.

Parameters:

asNumpy (bool=False) – If true, the values are returned as a numpy array instead of a list of Vec3s.

Returns:

The atomic positions

Return type:

list of Vec3 or np.ndarray

getBoxVectors(asNumpy: bool = False) List[Vec3] | ndarray

Get the periodic box vectors.

Parameters:

asNumpy (bool, optional, default=False) – If true, the values are returned as a numpy array instead of a list of Vec3s.

Returns:

The periodic box vectors.

Return type:

list of Vec3 or np.ndarray