5.20. AMBER PRMTOP topology parser

Reads an AMBER top file to build the system.

Amber keywords are turned into the following attributes:

AMBER flag MDAnalysis attribute
ATOM_NAME names
CHARGE charges
ATOMIC_NUMBER elements
MASS masses
BONDS_INC_HYDROGEN BONDS_WITHOUT_HYDROGEN bonds
ANGLES_INC_HYDROGEN ANGLES_WITHOUT_HYDROGEN angles
DIHEDRALS_INC_HYDROGEN DIHEDRALS_WITHOUT_HYDROGEN dihedrals / improper
ATOM_TYPE_INDEX type_indices
AMBER_ATOM_TYPE types
RESIDUE_LABEL resnames
RESIDUE_POINTER residues

Note

The Amber charge is converted to electron charges as used in MDAnalysis and other packages. To get back Amber charges, multiply by 18.2223.

Chamber-style Amber topologies (i.e. topologies generated via parmed conversion of a CHARMM topology to an AMBER one) are not currently supported. Support will likely be added in future MDAnalysis releases.

As of version 2.0.0, elements are no longer guessed if ATOMIC_NUMBER records are missing. In those scenarios, if elements are necessary, users will have to invoke the element guessers after parsing the topology file. Please see MDAnalysis.topology.guessers for more details.

5.20.1. Classes

class MDAnalysis.topology.TOPParser.TOPParser(filename)[source]

Reads topology information from an AMBER top file.

Reads the following Attributes if in topology: - Atomnames - Charges - Masses - Elements - Atomtypes - Resnames - Type_indices - Bonds - Angles - Dihedrals (inc. impropers)

The format is defined in PARM parameter/topology file specification. The reader tries to detect if it is a newer (AMBER 12?) file format by looking for the flag “ATOMIC_NUMBER”.

Notes

Elements are obtained from the atomic numbers (if present). If a given input atomic number does not belong to an element (usually either -1 or 0), the element will be assigned an empty record.

Changed in version 0.7.6: parses both amber10 and amber12 formats

Changed in version 0.19.0: parses bonds, angles, dihedrals, and impropers

Changed in version 1.0.0: warns users that chamber-style topologies are not currently supported

Changed in version 2.0.0: no longer guesses elements if missing

close()

Close the trajectory file.

convert_forces_from_native(force, inplace=True)

Conversion of forces array force from native to base units

Parameters:
  • force (array_like) – Forces to transform
  • inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data

Note

By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.

New in version 0.7.7.

convert_forces_to_native(force, inplace=True)

Conversion of force array force from base to native units.

Parameters:
  • force (array_like) – Forces to transform
  • inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data

Note

By default, the input force is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.

New in version 0.7.7.

convert_pos_from_native(x, inplace=True)

Conversion of coordinate array x from native units to base units.

Parameters:
  • x (array_like) – Positions to transform
  • inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data

Note

By default, the input x is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.

Changed in version 0.7.5: Keyword inplace can be set to False so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.

convert_pos_to_native(x, inplace=True)

Conversion of coordinate array x from base units to native units.

Parameters:
  • x (array_like) – Positions to transform
  • inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data

Note

By default, the input x is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.

Changed in version 0.7.5: Keyword inplace can be set to False so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.

convert_time_from_native(t, inplace=True)

Convert time t from native units to base units.

Parameters:
  • t (array_like) – Time values to transform
  • inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data

Note

By default, the input t is modified in place and also returned (although note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.) In-place operations improve performance because allocating new arrays is avoided.

Changed in version 0.7.5: Keyword inplace can be set to False so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.

convert_time_to_native(t, inplace=True)

Convert time t from base units to native units.

Parameters:
  • t (array_like) – Time values to transform
  • inplace (bool, optional) – Whether to modify the array inplace, overwriting previous data

Note

By default, the input t is modified in place and also returned. (Also note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.)

Changed in version 0.7.5: Keyword inplace can be set to False so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.

convert_velocities_from_native(v, inplace=True)

Conversion of velocities array v from native to base units

Parameters:
  • v (array_like) – Velocities to transform
  • inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data

Note

By default, the input v is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.

New in version 0.7.5.

convert_velocities_to_native(v, inplace=True)

Conversion of coordinate array v from base to native units

Parameters:
  • v (array_like) – Velocities to transform
  • inplace (bool (optional)) – Whether to modify the array inplace, overwriting previous data

Note

By default, the input v is modified in place and also returned. In-place operations improve performance because allocating new arrays is avoided.

New in version 0.7.5.

parse(**kwargs)[source]

Parse Amber PRMTOP topology file filename.

Returns:
Return type:A MDAnalysis Topology object
parse_bonded(num_per_record, numlines)[source]

Extracts bond information from PARM7 format files

Parameters:
  • num_per_record (int) – The number of entries for each record in section
  • numlines (int) – The number of lines to be parsed for this section

Note

For the bond/angle sections of parm7 files, the atom numbers are set to coordinate array index values. As detailed in http://ambermd.org/formats.html to recover the actual atom number, one should divide the values by 3 and add 1. Here, since we want to satisfy zero-indexing, we only divide by 3.

parse_charges(num_per_record, numlines)[source]

Extracts the partial charges for each atom

Parameters:
  • num_per_record (int) – The number of entries for each record in section (unused input)
  • numlines (int) – The number of lines to be parsed in current section
Returns:

attr – A Charges instance containing the partial charges of each atom as defined in the parm7 file

Return type:

Charges

parse_chunks(data, chunksize)[source]

Helper function to parse AMBER PRMTOP bonds/angles.

Parameters:
  • data (list of int) – Input list of the parm7 bond/angle section, zero-indexed
  • num_per_record (int) – The number of entries for each record in the input list
Returns:

vals – A list of tuples containing the atoms involved in a given bonded interaction

Return type:

list of int tuples

Note

In the parm7 format this information is structured in the following format: [ atoms 1:n, internal index ] Where 1:n represent the ids of the n atoms involved in the bond/angle and the internal index links to a given set of FF parameters. Therefore, to extract the required information, we split out the list into chunks of size num_per_record, and only extract the atom ids.

parse_dihedrals(diha, dihh)[source]

Combines hydrogen and non-hydrogen containing AMBER dihedral lists and extracts sublists for conventional dihedrals and improper angles

Parameters:
  • diha (list of tuples) – The atom ids of dihedrals not involving hydrogens
  • dihh (list of tuples) – The atom ids of dihedrals involving hydrogens
Returns:

  • dihedrals (Dihedrals) – A Dihedrals instance containing a list of all unique dihedrals as defined by the parm7 file
  • impropers (Impropers) – A Impropers instance containing a list of all unique improper dihedrals as defined by the parm7 file

Note

As detailed in http://ambermd.org/formats.html, the dihedral sections of parm7 files contain information about both conventional dihedrals and impropers. The following must be accounted for: 1) If the fourth atom in a dihedral entry is given a negative value, this indicates that it is an improper. 2) If the third atom in a dihedral entry is given a negative value, this indicates that it 1-4 NB interactions are ignored for this dihedrals. This could be due to the dihedral within a ring, or if it is part of a multi-term dihedral definition or if it is an improper.

parse_elements(num_per_record, numlines)[source]

Extracts the atomic numbers of each atom and converts to element type

Parameters:
  • num_per_record (int) – The number of entries for each record in section(unused input)
  • numlines (int) – The number of lines to be pasred in current section
Returns:

attr – A Elements instance containing the element of each atom as defined in the parm7 file

Return type:

Elements

Note

If the record contains unknown atomic numbers (e.g. <= 0), these will be treated as unknown elements and assigned an empty string value. See issues #2306 and #2651 for more details.

Changed in version 2.0.0: Unrecognised elements will now return a empty string. The parser will no longer attempt to guess the element by default.

parse_masses(num_per_record, numlines)[source]

Extracts the mass of each atom

Parameters:
  • num_per_record (int) – The number of entries for each record in section (unused input)
  • numlines (int) – The number of lines to be parsed in current section
Returns:

attr – A Masses instance containing the mass of each atom as defined in the parm7 file

Return type:

Masses

parse_names(num_per_record, numlines)[source]

Extracts atoms names from parm7 file

Parameters:
  • num_per_record (int) – The number of entries for each record in the section (unused input)
  • numlines (int) – The number of lines to be parsed in current section
Returns:

attr – A Atomnames instance containing the names of each atom as defined in the parm7 file

Return type:

Atomnames

parse_residx(num_per_record, numlines)[source]

Extracts the residue pointers for each atom

Parameters:
  • num_per_record (int) – The number of entries for each record in section (unused input)
  • numlines (int) – The number of lines to be parsed in current section
Returns:

vals – A list of zero-formatted residue pointers for each atom

Return type:

list of int

parse_resnames(num_per_record, numlines)[source]

Extracts the names of each residue

Parameters:
  • num_per_record (int) – The number of entries for each recrod in section (unused input)
  • numlines (int) – The number of lines to be parsed in current section
Returns:

attr – A Resnames instance containing the names of each residue as defined in the parm7 file

Return type:

Resnames

parse_type_indices(num_per_record, numlines)[source]

Extracts the index of atom types of the each atom involved in Lennard Jones (6-12) interactions.

Parameters:
  • num_per_record (int) – The number of entries for each record in section (unused input)
  • numlines (int) – The number of lines to be parsed in current section
Returns:

A TypeIndices instance containing the LJ 6-12 atom type index for each atom

Return type:

attr TypeIndices

parse_types(num_per_record, numlines)[source]

Extracts the force field atom types of each atom

Parameters:
  • num_per_record (int) – The number of entries for each record in section (unused input)
  • numlines (int) – The number of lines to be parsed in current section
Returns:

attr – A Atomtypes instance containing the atom types for each atom as defined in the parm7 file

Return type:

Atomtypes

parsesection_mapper(numlines, mapper)[source]

Parses FORTRAN formatted section, and returns a list of all entries in each line

Parameters:
  • numlines (int) – The number of lines to be parsed in this section
  • mapper (lambda operator) – Operator to format entries in current section
Returns:

section – A list of all entries in a given parm7 section

Return type:

list

skipper()[source]

TOPParser :class: helper function, skips lines of input parm7 file until we find the next %FLAG entry and return that

Returns:line – String containing the current line of the parm7 file
Return type:string