Bonds

The Chemical System contains a bonds property, which is an instance of the Bonds class. This contains the (possibly empty) bonding information on the system.

Creating and adding bonds

There are several ways of getting bonds into your chemical system. The following snippet will show some of the most common approaches:

from scm.libbase import ChemicalSystem
from scm.libbase import Bond

# Here we initialize a Chemical System from an XYZ file.
# The chemical system will not have any bonds at the moment, because
# bonds are not defined in the .xyz file.
water = ChemicalSystem.from_xyz("water.xyz")

# You can let the bond guessing algorithm simply guess the bonds:
water.guess_bonds()


# Or you can add them manually.
# First we clear the bonds we just guessed:
water.bonds.clear_bonds()

# and add bonds manually:
water.bonds.add_bond(0, 1, Bond(1.0)) # From the oxygen to the first Hydrogen
water.bonds.add_bond(0, 2, Bond(1.0)) # From the oxygen to the second Hydrogen


# You can also explicitly define the bonds in the System Block...
system_block = """
System
    Atoms
        O    0.0000000000000000  0.0000000000000000  0.3695041700000000
        H    0.0000000000000000  0.7838367199999999 -0.1847520900000000
        H    0.0000000000000000 -0.7838367199999999 -0.1847520900000000
    End
    BondOrders
        1 2 1.00
        1 3 1.00
    End
End
"""

# ...and then load the system:
water = ChemicalSystem(system_block)

# This is how you can iterate over the bonds and get the indices of the bonded atoms:
for i, j, bond in water.bonds:
    print(f"Atoms: {i}-{j} Bond order: {bond.order}")

Note on atom indexing: as you can also see from the example above, atom indexing from within python and in the System Block differ: In python the first atom has index 0, while in the System Block the first atom has index 1.

Bond Class

This contains the properties of a single bond.

class Bond
class Bond(order: float)
class Bond(order: float, lattice_displacements: ArrayLike)

A class representing a bond between two atoms.

get_inverted() Bond

Returns a copy of the bond with inverted lattice displacements.

property lattice_displacements: ndarray[Any, dtype[int64]]

Lattice displacements for the second atom that is part of the bond.

If atom A and B are bonded with lattice displacements [1, 0, 0], this means that atom A is bonded to the image of atom B that is displaced along the first lattice vector. Note that A > B in terms of atom indices.

property order: float

The order of the bond: 1 = single bond, 1.5 = aromatic bond, 2 = double bond, 3 = triple bond.

Bonds Class

class Bonds(num_atoms: int, num_lattice_vectors: int)

A class representing a set of bonds between atoms in a ChemicalSystem.

Methods to inspect the bonds:

num_bonds() int
num_bonds(atom: int) int
num_bonds(from_atom: int, to_atom: int) int

If called without an atom index, returns the total number of bonds between all atoms. If called with an atom index, returns the number of bonds of that atom. if called with two atom indices, returns the number of bonds between the two atoms.

get_bond_index(bond: Bond) int

“Returns the index of a bond in the .bonds attribute.

If the bond is not found (i.e. it is a bond contained in a different instance of Bonds), the number -1 is returned.

atoms_are_bonded(from_atom: int, to_atom: int) bool

Checks whether two atoms are bonded or not.

any_bond_between(from_atoms: List[int], to_atoms: List[int]) bool

Checks if there is any direct bond between two groups of atoms.

get_bonded_atoms(atom: int) List[int]

Returns a list of unique atom indices that are bonded to the given atom index.

get_atomic_valencies() ndarray[Any, dtype[float64]]

Returns the atomic valencies, i.e. the sum of all bond orders, for each atom.

property num_atoms: int

The total number of atoms in the system.

property num_lattice_vectors: int

The number of lattice vectors of the system.

Methods to iterate over bonds:

__len__() int

Returns the total number of bonds between all atoms.

__iter__() Iterator[Tuple[int, int, Bond]]

Returns an iterator over all bonds.

Unpacking the iterator yields a 3-tuple of two atom indices A and B, and a Bond instance. Note that A <= B to avoid double counting.

get_bonds_for_atom(from_atom: int) Iterator[Tuple[int, int, Bond]]

Returns an iterator over all bonds from a particular atom given by its index.

Unpacking the iterator yields a 3-tuple of two atom indices A and B, and a Bond instance. Note that A <= B.

get_bonds_between_atoms(from_atom: int, to_atom: int) Iterator[Tuple[int, int, Bond]]

Returns an iterator over all bonds from between two atoms given by their indices.

Unpacking the iterator yields a 3-tuple of two atom indices A and B, and a Bond instance. Note that A <= B.

Methods to modify the bonds:

add_bond(from_atom: int, to_atom: int, bond: Bond) None

Adds a bond between two atoms, given their indices.

add_bonds(from_atoms: ArrayLike, to_atoms: ArrayLike, bonds: List[Bond]) None
add_bonds(from_atoms: ArrayLike, to_atoms: ArrayLike, orders: ArrayLike) None

Adds bonds between multiple atom pairs at the same time. All arguments should have the same length.

remove_bond(bidx: int | Bond) None

Removes a bond, given either its index in the .bonds attribute or an instance of a Bond.

remove_bonds(bidx: List[int]) None

Removes multiple bonds at once, given their indices in the .bonds attribute

remove_bonds_to_atom(atom: int) None

Removes all bonds to a particular atom, given its index.

remove_bonds_to_atoms(remove_atoms: ArrayLike) None

Removes all bonds for a list of atoms given by their indices.

remove_bonds_between_atoms(from_atom: int, to_atom: int) None

Removes the bonds between two atoms, gived their indices.

clear_bonds() None

Removes all bonds.

Direct read-only access to the modified Compressed Sparse Row format:

property row_offset: ndarray[Any, dtype[int64]]

Bonds are stored as sparse matrix in a modified Compressed Spare Row format: The index in row_offset represents from_atom, and the value is an offset in column_index and bonds. The values in column_index represent to_atom. Because there can be duplicate indices (periodic bonds), it is not strictly a CSR matrix.

property column_index: ndarray[Any, dtype[int64]]

Bonds are stored as sparse matrix in a modified Compressed Spare Row format: The index in row_offset represents from_atom, and the value is an offset in column_index and bonds. The values in column_index represent to_atom. Because there can be duplicate indices (periodic bonds), it is not strictly a CSR matrix.

property bond_index: ndarray[Any, dtype[int64]]

Bonds are stored as sparse matrix in a modified Compressed Spare Row format: The index in row_offset represents from_atom, and the value is an offset in column_index and bonds. The values in column_index represent to_atom. Because there can be duplicate indices (periodic bonds), it is not strictly a CSR matrix.

property bonds: BondList

Bonds are stored as sparse matrix in a modified Compressed Spare Row format: The index in row_offset represents from_atom, and the value is an offset in column_index and bonds. The values in column_index represent to_atom. Because there can be duplicate indices (periodic bonds), it is not strictly a CSR matrix.

Supported file formats and conversions to other formats:

classmethod from_kf(kf: KFFile, section: str) Bonds

Constructs and returns new Bonds from a section on a KF file.

write_kf(kf: KFFile, section: str, write_list_format: bool = True) None

Writes Bonds to a section on a KF file.

The write_list_format argument determines whether the retrocompatible list format (with variables fromAtoms, toAtoms, bondOrders) is written. Otherwise the Compressed Sparse Row format is written directly.

classmethod from_list(num_atoms: int, num_lattice_vectors: int, from_atoms: ArrayLike, to_atoms: ArrayLike, bond_orders: ArrayLike, lattice_displacements: ArrayLike = [[]], indexing: int = 0) Bonds

Creates a set of bonds from lists of atom indices, bond orders and optional lattice displacements.

The indexing argument can be used to switch to one-based atom indices, as they are used in the AMS System block format.

classmethod from_sparse(num_lattice_vectors: int, row_offset: ArrayLike, column_index: ArrayLike, bond_orders: ArrayLike, lattice_displacements: ArrayLike = [[]], indexing: int = 0) Bonds

Low level construction of a set of bonds directly from the (modified) CRS storage.

The indexing argument can be used to switch to one-based atom indices, as they are used in the AMS System block format.

Bonds methods in the Chemical System

property ChemicalSystem.bonds: Bonds

The bonds of the system.

ChemicalSystem.has_bonds() bool

Checks if the system has any bonding information.

ChemicalSystem.clear_bonds() None

Removes all bonds from the system.

ChemicalSystem.guess_bonds() None

Guesses bonds based on the atomic elements and the geometry. Keeps existing bonds.

ChemicalSystem.set_lattice_displacements_from_minimum_image_convention() None

Applies the minimum image convention which sets the lattice_displacements for all bonds to be the shortest possible.

ChemicalSystem.num_molecules() int

Counts the number of connected molecules that are part of this system.

What is considered a molecule is based on the bonds of the system: anything connected by any bonds is considered to be a molecule. For a fully connected system, this method returns 1. If the system has no bonding information, each atom is its own molecule and this method returns the total number of atoms in the system.

ChemicalSystem.molecule_indices() ndarray[Any, dtype[int64]]

Returns a num_atoms sized array mapping the atoms to connected molecules of the system.

What is considered a molecule is based on the bonds of the system: anything connected by any bonds is considered to be a molecule. For a system consisting of two water molecules (where the first three atoms make up one of the molecules), this method returns [0, 0, 0, 1, 1, 1]. If the system has no bonding information, each atoms is considered its own molecule and this method returns [0, 1, 2, ..., num_atoms].

ChemicalSystem.split_into_molecules() List[ChemicalSystem]

Splits the system into individual molecules based on the connectivity between atoms.

The length of the returned list is num_molecules. Which set of atoms ends up in which element of the list is determined by the indices returned by the molecule_indices method.

The returned systems will have the same lattice as the original and a total charge of zero. Atomic attributes and regions are preserved, but the returned systems do not have any selected atoms.

This method is just a shortcut for: cs.split(cs.molecule_indices()).

ChemicalSystem.bond_cuts_molecule(from_atom: int, to_atom: int) bool
ChemicalSystem.bond_cuts_molecule(from_atom: Atom, to_atom: Atom) bool

Checks if removing the bonds between two atoms would cut the graph into two disjoint subgraphs.

ChemicalSystem.atom_is_in_ring(atom: int | Atom) bool

Checks if an atom is part of any ring.

An atom is in a ring, if removing the bond to (at least) one of its neighbors does not cut the graph into two disjoint subgraphs.

ChemicalSystem.molgraph_dijkstra(from_atidx: int, dist_func: molgraph_edgeweight_func_t, to_atidx: int = -1) Tuple[List[float], List[int], List[bool]]

General method implementing the Dijkstra algorithm on the molecular graph. The Dijktra algorithm solves the single-source shortest path problem in weighted graphs with non-negative weights.

  • from_atidx is the index of source atom, i.e. the atom from which we start searching for paths.

  • dist_func is a method (of type molgraph_edgeweight_func_t = Callable[[int, int, Bond], float]) taking two atomic indices and a Bond instance as input. It returns the weight of the graph edge, aka the distance between the two atoms. The returned distance may not be negative, but you can return float("inf") to ignore a particular edge in the search.

  • to_atidx marks an optional target atom and can be used to solve the single-pair shortest path problem by prematurely terminating the Dijkstra algorithm as soon as the shortest path to the target atom has been found.

Returns a 3-tuple consisting of:

  1. path_lengths is a list of length num_atoms that upon return of the method will be populated with the lengths of the shortest paths from the source atom. Unreachable atoms, aka atoms in other molecules, will have an infinite path length.

  2. path_through is a list of length num_atoms that upon return will for each atom contain the preceding atoms index in the shortest path from source atom. By stepping backwards through this array all shortest paths from the source atom can be recovered. If no path has been found to a particular atom, the corresponding list element will be set to -1.

  3. visited is a list of length num_atoms marking which atoms have been visited by the Dijkstra algorithm prior to its termination. An atom counts as visited if it has been popped of the top of the priority queue to have its neighbors added to the queue. Note that in case of premature termination of the algorithm due to using the to_atidx argument, shortest paths have only been determined for all visited atoms. (Other atoms may still have a path length set, but it is not guaranteed to be the shortest path from the source.)

ChemicalSystem.shortest_path_lengths_from(from_atom: int | Atom) List[float]
ChemicalSystem.shortest_path_lengths_from(from_atom: int | Atom, dist_func: molgraph_edgeweight_func_t) List[float]

Returns the lengths of all shortest paths from from_atom to all other atoms.

ChemicalSystem.shortest_path_between(from_atom: int, to_atom: int, dist_func: molgraph_edgeweight_func_t) List[int] | None
ChemicalSystem.shortest_path_between(from_atom: int, to_atom: int) List[int] | None
ChemicalSystem.shortest_path_between(from_atom: Atom, to_atom: Atom, dist_func: molgraph_edgeweight_func_t) List[int] | None
ChemicalSystem.shortest_path_between(from_atom: Atom, to_atom: Atom) List[int] | None

Returns the shortest path between two atoms as a list of atom indices.

Note that the list contains both from_atom and to_atom as the first and last element. The length of the path (measured in number of hops) is therefore equal to the length of the returned list minus one.

ChemicalSystem.shortest_path_length_between(from_atom: int, to_atom: int, dist_func: molgraph_edgeweight_func_t) float
ChemicalSystem.shortest_path_length_between(from_atom: int, to_atom: int) float
ChemicalSystem.shortest_path_length_between(from_atom: Atom, to_atom: Atom, dist_func: molgraph_edgeweight_func_t) float
ChemicalSystem.shortest_path_length_between(from_atom: Atom, to_atom: Atom) float

Returns the length of the shortest path between two atoms.