Atoms#

class autode.atoms.Atom(atomic_symbol: str, x: Any = 0.0, y: Any = 0.0, z: Any = 0.0, atom_class: int | None = None, partial_charge: float | None = None)#
__init__(atomic_symbol: str, x: Any = 0.0, y: Any = 0.0, z: Any = 0.0, atom_class: int | None = None, partial_charge: float | None = None)#

Atom class. Centered at the origin by default. Can be initialised from positional or keyword arguments:

>>> import autode as ade
>>> ade.Atom('H')
Atom(H, 0.0000, 0.0000, 0.0000)
>>>
>>> ade.Atom('H', x=1.0, y=1.0, z=1.0)
Atom(H, 1.0000, 1.0000, 1.0000)
>>>
>>> ade.Atom('H', 1.0, 1.0, 1.0)
Atom(H, 1.0000, 1.0000, 1.0000)
Parameters:
  • y – y coordinate in 3D space (Å)

  • z – z coordinate in 3D space (Å)

  • atom_class – Fictitious additional labels to distinguish otherwise identical atoms. Useful in finding bond isomorphisms over identity reactions

  • partial_charge – Partial atomic charge in units of e, determined by the atomic envrionment. Not an observable property.

property atomic_number: int#

Atomic numbers are the position in the elements (indexed from zero), plus one. Example:

>>> import autode as ade
>>> atom = ade.Atom('C')
>>> atom.atomic_number
6
property atomic_symbol: str#

A more interpretable alias for Atom.label. Should be present in the elements. Example:

>>> import autode as ade
>>> atom = ade.Atom('Zn')
>>> atom.atomic_symbol
'Zn'
property coord: Coordinate#

Position of this atom in space. Coordinate has attributes x, y, z for the Cartesian displacements. Example:

>>> import autode as ade
>>> atom = ade.Atom('H')
>>> atom.coord
Coordinate([0. 0. 0.] Å)

To initialise at a different position away from the origin

>>> ade.Atom('H', x=1.0).coord
Coordinate([1. 0. 0.] Å)
>>> ade.Atom('H', x=1.0).coord.x
1.0

Coordinates are instances of autode.values.ValueArray, so can be converted from the default angstrom units to e.g. Bohr

>>> ade.Atom('H', x=1.0, y=-1.0).coord.to('a0')
Coordinate([1.889  -1.889  0. ] bohr)
copy() Atom#
property covalent_radius: Distance#

Covalent radius for this atom. Example:

>>> import autode as ade
>>> ade.Atom('H').covalent_radius
Distance(0.31 Å)
property group: int#

Group of the periodic table is this atom in. 0 if not found. Example:

>>> import autode as ade
>>> ade.Atom('C').group
14
property is_metal: bool#

Is this atom a metal? Defines metals to be up to and including: Ga, Sn, Bi. Example:

>>> import autode as ade
>>> ade.Atom('C').is_metal
False
>>> ade.Atom('Zn').is_metal
True
is_pi(valency: int) bool#

Determine if this atom is a ‘π-atom’ i.e. is unsaturated. Only approximate! Example:

>>> import autode as ade
>>> ade.Atom('C').is_pi(valency=3)
True
>>> ade.Atom('H').is_pi(valency=1)
False
Return type:

(bool)

property mass: Mass#

Alias of weight. Returns Atom.weight an so can be converted to different units. For example, to convert the mass to electron masses:

>>> import autode as ade
>>> ade.Atom('H').mass.to('me')
Mass(1837.36222 m_e)
property maximal_valance: int#

The maximum/maximal valance that this atom supports in any charge state (most commonly). i.e. for H the maximal_valance=1. Useful for generating molecular graphs

property period: int#

Period of the periodic table is this atom in. 0 if not found. Example:

>>> import autode as ade
>>> ade.Atom('C').period
2
rotate(axis: ndarray | Sequence, theta: Angle | float, origin: ndarray | Sequence | None = None) None#

Rotate this atom theta radians around an axis given an origin. By default the rotation is applied around the origin with the angle in radians (unless an autode.values.Angle). Rotation is applied in place. To rotate a H atom around the z-axis:

>>> import autode as ade
>>> atom = ade.Atom('H', x=1.0)
>>> atom.rotate(axis=[0.0, 0.0, 1.0], theta=3.14)
>>> atom.coord
Coordinate([-1.  0.  0.] Å)

With an origin:

>>> import autode as ade
>>> atom = ade.Atom('H')
>>> atom.rotate(axis=[0.0, 0.0, 1.0], theta=3.14, origin=[1.0, 0.0, 0.0])
>>> atom.coord
Coordinate([2.  0.  0.] Å)

And with an angle not in radians:

>>> import autode as ade
>>> from autode.values import Angle
>>>
>>> atom = ade.Atom('H', x=1.0)
>>> atom.rotate(axis=[0.0, 0.0, 1.0], theta=Angle(180, units='deg'))
>>> atom.coord
Coordinate([-1.  0.  0.] Å)
Parameters:

origin – Rotate about this origin. shape = (3,) if no origin is specified then the atom is rotated without translation.

property tm_row: int | None#

Row of transition metals that this element is in. Returns None if this atom is not a metal. Example:

>>> import autode as ade
>>> ade.Atom('Zn').tm_row
1
translate(*args, **kwargs) None#

Translate this atom by a vector in place. Arguments should be coercible into a coordinate (i.e. length 3). Example:

>>> import autode as ade
>>> atom = ade.Atom('H')
>>> atom.translate(1.0, 0.0, 0.0)
>>> atom.coord
Coordinate([1. 0. 0.] Å)

Atoms can also be translated using numpy arrays:

>>> import autode as ade
>>> import numpy as np
>>>
>>> atom = ade.Atom('H')
>>> atom.translate(np.ones(3))
>>> atom.coord
Coordinate([1. 1. 1.] Å)
>>>
>>> atom.translate(vec=-atom.coord)
>>> atom.coord
Coordinate([0. 0. 0.] Å)
Keyword Arguments:

vec (np.ndarray) – Shape = (3,)

property vdw_radius: Distance#

Van der Waals radius for this atom. Example:

>>> import autode as ade
>>> ade.Atom('H').vdw_radius
Distance(1.1 Å)
property weight: Mass#

Atomic weight. Example:

>>> import autode as ade
>>> ade.Atom('C').weight
Mass(12.0107 amu)
>>>
>>> ade.Atom('C').weight == ade.Atom('C').mass
True
class autode.atoms.AtomCollection(atoms: List[Atom] | Atoms | None = None)#
__init__(atoms: List[Atom] | Atoms | None = None)#

Collection of atoms, used as a base class for a species, complex or transition state.

angle(i: int, j: int, k: int) Angle#

Angle between three atoms i-j-k, where the atoms are indexed from zero:

E_i  --- E_j
            \
      θ      E_k

Example:

>>> from autode import Atom, Molecule
>>> h2o = Molecule(atoms=[Atom('H', x=-1), Atom('O'), Atom('H', x=1)])
>>> h2o.angle(0, 1, 2).to('deg')
Angle(180.0 °)
Parameters:
  • j (int) – — middle

  • k (int) – — right

Returns:

Angle

Return type:

(autode.values.Angle)

Raises:

(ValueError) – If any of the atom indexes are not present

property atoms: Atoms | None#

Constituent atoms of this collection

property com: Coordinate | None#

Centre of mass of this atom collection

Raises:

(ValueError) – If there are no atoms

property coordinates: Coordinates | None#

Numpy array of coordinates

dihedral(w: int, x: int, y: int, z: int) Angle#

Dihedral angle between four atoms (x, y, z, w), where the atoms are indexed from zero:

E_w  --- E_x
            \  φ
             \
              E_y ---- E_z

Example:

>>> from autode import Atom, Molecule
>>> h2s2 = Molecule(atoms=[Atom('S',  0.1527, 0.9668, -0.9288),
...                        Atom('S',  2.0024, 0.0443, -0.4227),
...                        Atom('H', -0.5802, 0.0234, -0.1850),
...                        Atom('H',  2.1446, 0.8424,  0.7276)])
>>> h2s2.dihedral(2, 0, 1, 3).to('deg')
Angle(-90.0 °)
Parameters:
  • x (int) – – second –

  • y (int) – – third –

  • z (int) – – fourth –

Returns:

Dihedral angle

Return type:

(autode.values.Angle)

Raises:

(ValueError) – If any of the atom indexes are not present in the molecule

distance(i: int, j: int) Distance#
eqm_bond_distance(i: int, j: int) Distance#
property mass: Mass#

Molecular weight

property moi: MomentOfInertia | None#

Moment of inertia matrix (I)

property n_atoms: int#

Number of atoms in this collection

property weight: Mass#

Molecular weight

class autode.atoms.Atoms(iterable=(), /)#
are_linear(angle_tol: ~autode.values.Angle = Angle(1.0 °)) bool#

Are these set of atoms colinear?

Returns:

Whether the atoms are linear

Return type:

(bool)

are_planar(distance_tol: ~autode.values.Distance = Distance(0.001 Å)) bool#

Do all the atoms in this set lie in a single plane?

Return type:

(bool)

property com: Coordinate#

Centre of mass of these coordinates

\[\text{COM} = \frac{1}{M} \sum_i m_i R_i\]

where M is the total mass, m_i the mass of atom i and R_i it’s coordinate

property contain_metals: bool#

Do these atoms contain at least a single metal atom?

property coordinates: Coordinates#
copy() Atoms#

Copy these atoms, deeply

distance(i: int, j: int) Distance#

Distance between two atoms (Å), indexed from 0.

>>> import autode as ade
>>> mol = ade.Molecule(atoms=[ade.Atom('H'), ade.Atom('H', x=1.0)])
>>> mol.distance(0, 1)
Distance(1.0 Å)
Parameters:

j (int) – Atom index of the second atom

Returns:

Distance

Return type:

(autode.values.Distance)

Raises:

(ValueError)

eqm_bond_distance(i: int, j: int) Distance#

Equilibrium distance between two atoms. If known then use the experimental dimer distance, otherwise estimate if from the covalent radii of the two atoms. Example

Example:

>>> import autode as ade
>>> mol = ade.Molecule(atoms=[ade.Atom('H'), ade.Atom('H')])
>>> mol.distance(0, 1)
Distance(0.0 Å)
>>> mol.eqm_bond_distance(0, 1)
Distance(0.741 Å)
idxs_are_present(*args: int) bool#

Are all these indexes present in this set of atoms

property moi: MomentOfInertia#

Moment of inertia matrix (I):

    (I_00   I_01   I_02)
I = (I_10   I_11   I_12)
    (I_20   I_21   I_22)
Return type:

(autode.values.MomentOfInertia)

nvector(i: int, j: int) ndarray#

Normalised vector from atom i to atom j

Parameters:

j (int) –

Return type:

(np.ndarray)

Raises:

(IndexError) – If i or j are not present

remove_dummy() None#

Remove all the dummy atoms from this list of atoms

vector(i: int, j: int) ndarray#

Vector from atom i to atom j

Parameters:

j (int) –

Return type:

(np.ndarray)

Raises:

(IndexError) – If i or j are not present

class autode.atoms.DummyAtom(x, y, z)#
__init__(x, y, z)#

Dummy atom

Parameters:
  • y (float) – y

  • z (float) – z

property atomic_number#

The atomic number is defined as 0 for a dummy atom

property covalent_radius: Distance#

Dummy atoms have no radius

property mass: Mass#

Dummy atoms do not have any weight/mass

property vdw_radius: Distance#

Dummy atoms have no radius

property weight: Mass#

Dummy atoms do not have any weight/mass

class autode.atoms.PeriodicTable#
actinides = array(['Ac', 'Th', 'Pa', 'U', 'Np', 'Pu', 'Am', 'Cm', 'Bk', 'Cf', 'Es',        'Fm', 'Md', 'No', 'Lr'], dtype='<U2')#
actinoids = array(['Ac', 'Th', 'Pa', 'U', 'Np', 'Pu', 'Am', 'Cm', 'Bk', 'Cf', 'Es',        'Fm', 'Md', 'No', 'Lr'], dtype='<U2')#
classmethod element(period: int, group: int)#

Element given it’s index in the periodic table, excluding lanthanides and actinides.

Parameters:

group (int) –

Returns:

Atomic symbol of the element

Return type:

(str)

Raises:

(IndexError) – If such an element does not exist

classmethod group(n: int)#

Group of the periodic table, with 1 being the first period

Return type:

(np.ndarray(str))

Raises:

(ValueError) – If n is not valid group index

lanthanides = array(['La', 'Ce', 'Pr', 'Nd', 'Pm', 'Sm', 'Eu', 'Gd', 'Tb', 'Dy', 'Ho',        'Er', 'Tm', 'Yb', 'Lu'], dtype='<U2')#
lanthanoids = array(['La', 'Ce', 'Pr', 'Nd', 'Pm', 'Sm', 'Eu', 'Gd', 'Tb', 'Dy', 'Ho',        'Er', 'Tm', 'Yb', 'Lu'], dtype='<U2')#
classmethod period(n: int)#

Period of the periodic table, with 1 being the first period

Return type:

(np.ndarray(str))

Raises:

(ValueError) – If n is not valid period index

table = array([['H', '', '', '', '', '', '', '', '', '', '', '', '', '', '', '',         '', 'He'],        ['Li', 'Be', '', '', '', '', '', '', '', '', '', '', 'B', 'C',         'N', 'O', 'F', 'Ne'],        ['Na', 'Mg', '', '', '', '', '', '', '', '', '', '', 'Al', 'Si',         'P', 'S', 'Cl', 'Ar'],        ['K', 'Ca', 'Sc', 'Ti', 'V', 'Cr', 'Mn', 'Fe', 'Co', 'Ni', 'Cu',         'Zn', 'Ga', 'Ge', 'As', 'Se', 'Br', 'Kr'],        ['Rb', 'Sr', 'Y', 'Zr', 'Nb', 'Mo', 'Tc', 'Ru', 'Rh', 'Pd', 'Ag',         'Cd', 'In', 'Sn', 'Sb', 'Te', 'I', 'Xe'],        ['Cs', 'Ba', '', 'Hf', 'Ta', 'W', 'Re', 'Os', 'Ir', 'Pt', 'Au',         'Hg', 'Tl', 'Pb', 'Bi', 'Po', 'At', 'Rn'],        ['Fr', 'Ra', '', 'Rf', 'Db', 'Sg', 'Bh', 'Hs', 'Mt', 'Ds', 'Rg',         'Cn', 'Nh', 'Fl', 'Mc', 'Lv', 'Ts', 'Og']], dtype='<U2')#
classmethod transition_metals(row: int)#

Collection of transition metals (TMs) of a defined row. e.g.

row = 1 -> [Sc, Ti .. Zn]

Return type:

(np.ndarray(str))

Raises:

(ValueError) – If the row is not valid

autode.atoms.vdw_radii = {'Ac': 2.47, 'Ag': 2.11, 'Al': 1.84, 'Am': 2.44, 'Ar': 1.88, 'As': 1.85, 'At': 2.02, 'Au': 2.14, 'B': 1.92, 'Ba': 2.68, 'Be': 1.53, 'Bi': 2.07, 'Bk': 2.44, 'Br': 1.85, 'C': 1.7, 'Ca': 2.31, 'Cd': 2.18, 'Ce': 2.42, 'Cf': 2.45, 'Cl': 1.75, 'Cm': 2.45, 'Co': 2.0, 'Cr': 2.06, 'Cs': 3.43, 'Cu': 1.96, 'Dy': 2.31, 'Er': 2.29, 'Es': 2.45, 'Eu': 2.35, 'F': 1.47, 'Fe': 2.04, 'Fm': 2.45, 'Fr': 3.48, 'Ga': 1.87, 'Gd': 2.34, 'Ge': 2.11, 'H': 1.1, 'He': 1.4, 'Hf': 2.23, 'Hg': 2.23, 'Ho': 2.3, 'I': 1.98, 'In': 1.93, 'Ir': 2.13, 'K': 2.75, 'Kr': 2.02, 'La': 2.43, 'Li': 1.82, 'Lr': 2.46, 'Lu': 2.24, 'Md': 2.46, 'Mg': 1.73, 'Mn': 2.05, 'Mo': 2.17, 'N': 1.55, 'Na': 2.27, 'Nb': 2.18, 'Nd': 2.39, 'Ne': 1.54, 'Ni': 1.97, 'No': 2.46, 'Np': 2.39, 'O': 1.52, 'Os': 2.16, 'P': 1.8, 'Pa': 2.43, 'Pb': 2.02, 'Pd': 2.1, 'Pm': 2.38, 'Po': 1.97, 'Pr': 2.4, 'Pt': 2.13, 'Pu': 2.43, 'Ra': 2.83, 'Rb': 3.03, 'Re': 2.16, 'Rh': 2.1, 'Rn': 2.2, 'Ru': 2.13, 'S': 1.8, 'Sb': 2.06, 'Sc': 2.15, 'Se': 1.9, 'Si': 2.1, 'Sm': 2.36, 'Sn': 2.17, 'Sr': 2.49, 'Ta': 2.22, 'Tb': 2.33, 'Tc': 2.16, 'Te': 2.06, 'Th': 2.45, 'Ti': 2.11, 'Tl': 1.96, 'Tm': 2.27, 'U': 2.41, 'V': 2.07, 'W': 2.18, 'Xe': 2.16, 'Y': 2.32, 'Yb': 2.26, 'Zn': 2.01, 'Zr': 2.23}#

Although a π-bond may not be well defined, it is useful to have a notion of a bond about which there is restricted rotation. The below sets are used to define which atoms may be π-bonded to another