TOPContainer

A class for reading and GROMACS .top topology files.

Index

TOPContainer(*[, defaults, atomtypes, ...])

A class for managing GROMACS .top topology files.

TOPContainer.from_file(path)

Construct a new TOPContainer from the passed file path.

TOPContainer.to_file(path)

Construct a new .top file from this instance.

TOPContainer.allclose(other, *[, rtol, ...])

Return whether two TOPContainers are equivalent within a given tolerance.

TOPContainer.generate_pairs([func])

Construct and populate the pairs directive with explicit 1,4-pairs based on the available bonds.

TOPContainer.generate_pairs_nb([func])

Construct and populate the pairs_nb directive with explicit nonbonded pairs based on the available non-bonded atoms.

TOPContainer.copy([deep])

Return a copy of this instance.

TOPContainer.concatenate

Namespace with functions for adding new directive-specific rows.

API

class FOX.TOPContainer(*, defaults=None, atomtypes=None, moleculetype=None, atoms=None, system=None, molecules=None, bondtypes=None, pairtypes=None, angletypes=None, dihedraltypes=None, constrainttypes=None, nonbond_params=None, pairs=None, pairs_nb=None, bonds=None, angles=None, dihedrals=None)[source]

A class for managing GROMACS .top topology files.

Examples

>>> from FOX import TOPContainer

>>> input_file: str = ...
>>> output_file: str = ...

>>> top = TOPContainer.from_file(input_file)
>>> top.to_file(output_file)
DF_DTYPES = mappingproxy({'defaults': dtype([('nbfunc', '<i8'), ('comb_rule', '<i8'), ('gen_pairs', '<U3'), ('fudgeLJ', '<f8'), ('fudgeQQ', '<f8')]), 'atomtypes': dtype([('atom_type', '<U5'), ('atnum', '<i8'), ('mass', '<f8'), ('charge', '<f8'), ('particle_type', '<U1'), ('sigma', '<f8'), ('epsilon', '<f8')]), 'moleculetype': dtype([('molecule', 'O'), ('n_rexcl', '<i8')]), 'atoms': dtype([('molecule', 'O'), ('atom1', '<i8'), ('atom_type', '<U5'), ('res_num', '<i8'), ('res_name', '<U5'), ('atom_name', '<U5'), ('charge_group', '<i8'), ('charge', '<f8'), ('mass', '<f8')]), 'system': dtype([('name', 'O')]), 'molecules': dtype([('molecule', 'O'), ('n_mol', '<i8')]), 'bonds': dtype([('molecule', 'O'), ('atom1', '<i8'), ('atom2', '<i8'), ('func', '<i8')]), 'angles': dtype([('molecule', 'O'), ('atom1', '<i8'), ('atom2', '<i8'), ('atom3', '<i8'), ('func', '<i8')]), 'dihedrals': dtype([('molecule', 'O'), ('atom1', '<i8'), ('atom2', '<i8'), ('atom3', '<i8'), ('atom4', '<i8'), ('func', '<i8')]), 'pairs': dtype([('molecule', 'O'), ('atom1', '<i8'), ('atom2', '<i8'), ('func', '<i8')]), 'pairs_nb': dtype([('molecule', 'O'), ('atom1', '<i8'), ('atom2', '<i8'), ('func', '<i8')])})

A mapping holding the data types of all mandatory directives.

DF_DICT_DTYPES = mappingproxy({'pairtypes': mappingproxy({1: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('sigma', '<f8'), ('epsilon', '<f8')]), 2: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('fudgeQQ', '<f8'), ('qi', '<f8'), ('qj', '<f8'), ('sigma', '<f8'), ('epsilon', '<f8')])}), 'bondtypes': mappingproxy({1: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('b0', '<f8'), ('k', '<f8')]), 2: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('b0', '<f8'), ('k', '<f8')]), 3: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('b0', '<f8'), ('D', '<f8'), ('beta', '<f8')]), 4: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('b0', '<f8'), ('C', '<f8')]), 5: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8')]), 6: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('b0', '<f8'), ('k', '<f8')]), 7: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('b0', '<f8'), ('k', '<f8')]), 8: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('table_num', '<i8'), ('k', '<f8')]), 9: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('table_num', '<i8'), ('k', '<f8')]), 10: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('low', '<f8'), ('up', '<f8'), ('k', '<f8')])}), 'angletypes': mappingproxy({1: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('func', '<i8'), ('theta', '<f8'), ('k', '<f8')]), 2: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('func', '<i8'), ('theta', '<f8'), ('k', '<f8')]), 3: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('func', '<i8'), ('r1', '<f8'), ('r2', '<f8'), ('k', '<f8')]), 4: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('func', '<i8'), ('r1', '<f8'), ('r2', '<f8'), ('r3', '<f8'), ('k', '<f8')]), 5: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('func', '<i8'), ('theta', '<f8'), ('ktheta', '<f8'), ('ub0', '<f8'), ('kub', '<f8')]), 6: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('func', '<i8'), ('theta', '<f8'), ('C', '<f8')]), 8: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('func', '<i8'), ('table_num', '<i8'), ('k', '<f8')]), 9: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('func', '<i8'), ('a', '<f8'), ('k', '<f8')]), 10: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('func', '<i8'), ('theta', '<f8'), ('k', '<f8')])}), 'dihedraltypes': mappingproxy({1: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('atom4', '<U5'), ('func', '<i8'), ('phi0', '<f8'), ('k', '<f8'), ('n', '<i8')]), 2: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('atom4', '<U5'), ('func', '<i8'), ('xi0', '<f8'), ('k', '<f8')]), 3: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('atom4', '<U5'), ('func', '<i8'), ('C0', '<f8'), ('C1', '<f8'), ('C2', '<f8'), ('C3', '<f8'), ('C4', '<f8'), ('C5', '<f8')]), 4: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('atom4', '<U5'), ('func', '<i8'), ('phi0', '<f8'), ('k', '<f8'), ('n', '<i8')]), 5: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('atom4', '<U5'), ('func', '<i8'), ('C1', '<f8'), ('C2', '<f8'), ('C3', '<f8'), ('C4', '<f8'), ('C5', '<f8')]), 8: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('atom4', '<U5'), ('func', '<i8'), ('table_num', '<i8'), ('k', '<f8')]), 9: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('atom4', '<U5'), ('func', '<i8'), ('phi0', '<f8'), ('k', '<f8'), ('n', '<i8')]), 10: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('atom4', '<U5'), ('func', '<i8'), ('phi0', '<f8'), ('k', '<f8')]), 11: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('atom3', '<U5'), ('atom4', '<U5'), ('func', '<i8'), ('k', '<f8'), ('a0', '<f8'), ('a1', '<f8'), ('a2', '<f8'), ('a3', '<f8')])}), 'constrainttypes': mappingproxy({1: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8')]), 2: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8')])}), 'nonbond_params': mappingproxy({1: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('sigma', '<f8'), ('epsilon', '<f8')]), 2: dtype([('atom1', '<U5'), ('atom2', '<U5'), ('func', '<i8'), ('a', '<f8'), ('b', '<f8'), ('c', '<f8')])})})

A mapping holding the data types of all optional (dictionary of dataframe based) directives.

defaults

A dataframe holding the defaults directive.

atomtypes

A dataframe holding the atomtypes directive.

bondtypes

A dictionary of dataframes holding the bondtypes directive.

pairtypes

A dictionary of dataframes holding the pairtypes directive.

angletypes

A dictionary of dataframes holding the angletypes directive.

dihedraltypes

A dictionary of dataframes holding the dihedraltypes directive.

constrainttypes

A dictionary of dataframes holding the constrainttypes directive.

nonbond_params

A dictionary of dataframes holding the nonbond_params directive.

moleculetype

A dataframe holding the moleculetype directive.

atoms

A dataframe holding the atoms directive.

pairs

A dataframe holding the pairs directive.

bonds

A dataframe holding the bonds directive.

angles

A dataframe holding the angles directive.

dihedrals

A dataframe holding the dihedrals directive.

system

A dataframe holding the system directive.

molecules

A dataframe holding the molecules directive.

classmethod TOPContainer.from_file(path)[source]

Construct a new TOPContainer from the passed file path.

Parameters:

path (path-like object) – The path to the .top file

Returns:

A newly constructed .top container

Return type:

FOX.TOPContainer

TOPContainer.to_file(path)[source]

Construct a new .top file from this instance.

Parameters:

path (path-like object) – The path of the to-be created .top file

TOPContainer.allclose(other, *, rtol=1e-05, atol=1e-08, equal_nan=True)[source]

Return whether two TOPContainers are equivalent within a given tolerance.

Parameters:
  • other (TOPContainer) – The to-be compared TOPContainer

  • rtol (float) – The relative tolerance parameter (see Notes).

  • atol (float) – The absolute tolerance parameter (see Notes).

  • equal_nan (bool) – Whether to compare NaN’s as equal. If True, NaN’s in a will be considered equal to NaN’s in b in the output array.

Returns:

Whether the two containers are equivalent within a given tolerance.

Return type:

bool

See also

numpy.allclose()

Returns True if two arrays are element-wise equal within a tolerance.

TOPContainer.generate_pairs(func=1)[source]

Construct and populate the pairs directive with explicit 1,4-pairs based on the available bonds.

Parameters:

func ({1, 2}) – The func type as used for the new pairs.

TOPContainer.generate_pairs_nb(func=1)[source]

Construct and populate the pairs_nb directive with explicit nonbonded pairs based on the available non-bonded atoms.

Parameters:

func ({1}) – The func type as used for the new pairs.

TOPContainer.copy(deep=True)[source]

Return a copy of this instance.

Parameters:

deep (bool) – Whether a deep copy should be created or not

Return type:

A copy of this instance

TOPContainer.concatenate

Namespace with functions for adding new directive-specific rows.

atomtypes(self, *, atnum=None, symbol=None, atom_type=None, charge=0.0, sigma=0.0, epsilon=0.0, particle_type='A')

Add one or more atom types to the atomtypes directive.

Examples

>>> import FOX

>>> top: FOX.TOPContainer = ...
>>> top.concatenate.atomtypes(atnum=[6, 6, 7], sigma: [1.5, 1.2, 5.0], charge=0)
Parameters:
  • atnum/symbol (array-like) – One or more atomic numbers _or_ atomic symbols

  • atom_type (array-like) – One or more atom types. If not provided use normal atomic symbols instead

  • charge (array-like) – One or more atomic charges

  • sigma (array-like) – One or more Lennard-Jones sigma values

  • epsilon (array-like) – One or more Lennard-Jones epsilon values

  • particle_type ({"A", "S", "V", "D"}) – One or more particule types

FOX.io._top_concat._TOPConcat.nonbond_params(self, atom1, atom2, *, func, **kwargs)

Add one or more atom types to the nonbond_params directive.

Examples

>>> import FOX

>>> top: FOX.TOPContainer = ...
>>> top.concatenate.nonbond_params(
...     ["C12", "C12"], ["C12", "H33"], func=1, epsilon=1.0, sigma=0.5
... )
Parameters:
  • atom1 (array-like) – One or more atom types for the first atom defining the bond

  • atom2 (array-like) – One or more atom types for the second atom defining the bond

  • func ({1, 2}) – The type of potential function for the non-bonded potential, 1 representing Lennard-Jones and 2 Buckingham

  • **kwargs (array-like) – Func-specific extra (optional) arguments: * 1: sigma and epsilon * 2: a, b and c

FOX.io._top_concat._TOPConcat.atoms(self, atom_type, molecule, *, res_num, res_name, atom1=None, atom_name=None, charge_group=None)

Add one or more atom types to the atoms directive.

Examples

>>> import FOX

>>> top: FOX.TOPContainer = ...
>>> top.concatenate.atoms(
...     molecule="mol1", res_num=5, res_name="OLA",
...     atom_type=["C12", "C12", "C13", "O38"],
... )
Parameters:
  • molecule (array-like) – One or more molecule names; must be present in the moleculetype directive

  • res_num (array-like) – One or more residue numbers

  • res_name (array-like) – One or more residue names

  • atom_type (array-like) – One or more atom types

  • atom1 (Array-like) – One or more atomic indices for the new atoms. Automatically inferred if unspecified.

  • atom_name (array-like) – One or more atom names. Defaults to the same values as atom_type if unspecified.

  • charge_group (array-like) – One or more charge groups. Defaults to the atomic index of unspecified.

FOX.io._top_concat._TOPConcat.pairs(self, atom1, atom2, molecule, *, func)

Add one or more atom types to the pairs directive.

Examples

>>> import FOX

>>> top: FOX.TOPContainer = ...
>>> top.concatenate.pairs([1, 2, 3], [2, 3, 1], func=1)
Parameters:
  • molecule (array-like) – One or more molecule names; must be present in the moleculetype directive

  • atom1 (array-like) – One or more atomic indices for the first atom defining the bond

  • atom2 (array-like) – One or more atomic indices for the second atom defining the bond

  • func ({1, 2}) – The type of potential function for the non-bonded potential, 1 representing Lennard-Jones and 2 Buckingham

FOX.io._top_concat._TOPConcat.pairs_nb(self, atom1, atom2, molecule, *, func=1)

Add one or more atom types to the pairs_nb directive.

Examples

>>> import FOX

>>> top: FOX.TOPContainer = ...
>>> top.concatenate.pairs_nb([1, 2, 3], [2, 3, 1], func=1)
Parameters:
  • molecule (array-like) – One or more molecule names; must be present in the moleculetype directive

  • atom1 (array-like) – One or more atomic indices for the first atom defining the bond

  • atom2 (array-like) – One or more atomic indices for the second atom defining the bond

  • func ({1}) – The type of potential function for the non-bonded potential, 1 representing Lennard-Jones