API

Reaction Templates

scine_art.reaction_template.MoFrAtIndices: TypeAlias = Tuple[int, int, int]

A triplet of integers pointing at a specific atom within a template. Templates group atoms into fragments and fragments into molecule, this triplet encodes this grouping (<mol_idx>, <frag_idx>, <atom_idx>). The first integer indicating the index of the molecule the atom is in, the second index pointing to the fragment within the molecule, the third pointing at an atom within this fragment.

scine_art.reaction_template.MoFrAtPair: TypeAlias = Tuple[MoFrAtIndices, MoFrAtIndices]

Matches two atoms. Both atoms are encoded based on the molecule-fragment-atom grouping.

scine_art.reaction_template.BondChangesContainer: TypeAlias = SideContainer[AssosDissos[MoFrAtPair]]

A nested dictionary encoding all bond associations ('assos') and dissociations ('dissos') that are part of a reaction template. Bond modifications are given per side of the reaction (outer dictionary key).

Example

>>> bond_changes: BondChangesContainer = ...
>>> (mol_idx1, frag_idx1, atom_idx1), (mol_idx2, frag_idx2, atom_idx2) = \
>>>      bond_changes['lhs']['assos']
scine_art.reaction_template.SideConversionTree: TypeAlias = List[List[List[MoFrAtIndices]]]

Maps one atom encoded by a molecule-fragment-atom index set onto another that is encoded the same way.

Example

>>> mapping: SideConversionTree = ...
>>> rhs_mol_idx, rhs_frag_idx, rhs_atom_idx = \
>>>      mapping[lhs_mol_idx][lhs_frag_idx][lhs_atom_idx]
scine_art.reaction_template.ShapesByMoFrAt: TypeAlias = List[List[List[Optional[masm.shapes.Shape]]]]

Lists atom shapes for atoms encoded by a molecule-fragment-atom triple

Example

>>> shapes: ShapesByMoFrAt = ...
>>> atom_shape = shapes[mol_idx][frag_idx][atom_idx]
scine_art.reaction_template.Match: TypeAlias = Dict[str, List[Tuple[Tuple[int, int], Tuple[int, int]]]]

Indices references the given lists of molecules in the generating call.

Example

>>> molecules = [masm_mol1, masm_mol2]
>>> match: Match = matching_function(molecules, ...)
>>> for association in match['assos']:
>>>     mol_idx1, atom_idx1 = association[0]
>>>     mol_idx2, atom_idx2 = association[1]
>>>     # connect the specified atoms to build product structures
>>>     #   (repeat with disconnections for dissociations)
class scine_art.reaction_template.AssosDissos(assos, dissos)[source]

A dictionary like container with the keys assos and dissos indicating data based on forming and breaking bonds.

class scine_art.reaction_template.ReactionTemplate(fragments, lhs_rhs_mappings, shapes, barriers=None, elementary_step_id=None)[source]

A class encoding a single reaction template, with a basic constructor, requiring detailed data structures.

For simpler methods of constructing a reaction template, please see the static methods of this class: from_trajectory_spline, from_trajectory and from_aligned_structures.

Parameters:
fragmentsTemplateTypeContainer[List[List[masm.Molecule]]]

The fragments on both sides of the reaction. Given once for each template flavor/type to be held.

lhs_rhs_mappingsTemplateTypeContainer[SideConversionTree]

The mapping of the atom between lhs and rhs of the templates.

shapesTemplateTypeContainer[ShapesByMoFrAt]

The shapes of the ligand field at all atoms in given molecular fragments.

barriersOptional[Tuple[float, float]], optional

The lowest known barriers in kJ/mol for forward. ('lhs' -> 'rhs') and backward ('rhs' -> 'lhs') reaction.

elementary_step_idOptional[str], optional

The ID of the/an analyzed elementary step in a SCINE Database.

apply(molecules, energy_cutoff=300, enforce_atom_shapes=True, enforce_forward=False, enforce_backward=False, allowed_directions=None)[source]

Tries to apply a given template to the presented molecules.

Parameters:
moleculesList[masm.Molecule]

The molecules to apply the template to.

energy_cutofffloat, optional

Only apply the template if a barrier in the trialed direction of less than this value has been reported, by default 300 (kJ/mol)

enforce_atom_shapesbool, optional

If true only allow atoms with the same coordination sphere shapes to be considered matching, by default True.

enforce_forwardbool, optional

Enforce the trial of the forward direction independent of the energy_cutoff, by default False

enforce_backwardbool, optional

Enforce the trial of the backward direction independent of the energy_cutoff, by default False

allowed_directionsList[str], optional

Filter the allowed directions to be trial, by default None, indicating both forward and backward directions are allowed.

Returns:
List[Tuple[List[masm.Molecule], List[List[Tuple[int, int]]], masm.Molecule, List[List[Tuple[int, int]]], Match]]

Results of all possible applications of the given template, one tuple for each possibility. Each tuple contains: 0. The list of resulting molecules 1. The mapping of all atoms from given molecules to those in the resulting products (0.). 2. The fused transition state, where all associative bond modifications are applied but non of the dissociate ones. 3. The mapping of all atoms from given molecules to those in the fused transition state (2.).

rtype:

List[Tuple[List[Molecule], List[List[Tuple[int, int]]], Molecule, List[List[Tuple[int, int]]], Dict[str, List[Tuple[Tuple[int, int], Tuple[int, int]]]]]] ..

determine_all_matches(molecules, energy_cutoff=300, enforce_atom_shapes=True, enforce_lhs=False, enforce_rhs=False, allowed_sides=None)[source]

Tries to find matching fragment patterns of the reaction template in the given molecules.

Parameters:
moleculesList[masm.Molecule]

The molecules to apply the template to.

energy_cutofffloat, optional

Only apply the template if a barrier in the trialed direction of less than this value has been reported, by default 300 (kJ/mol)

enforce_atom_shapesbool, optional

If true only allow atoms with the same coordination sphere shapes to be considered matching, by default True.

enforce_lhsbool, optional

Enforce the trial of the left-hand side fragments of the reaction template independent of the energy_cutoff, by default False

enforce_rhsbool, optional

Enforce the trial of the right-hand side fragments of the reaction template independent of the energy_cutoff, by default False

allowed_sidesOptional[List[str]], optional

Filter the allowed sides to be tested, by default None, indicating both left-hand side (lhs) and right-hand side (rhs) are allowed.

Returns:
Optional[List[Match]]

A list of matches that were found.

rtype:

Optional[List[Dict[str, List[Tuple[Tuple[int, int], Tuple[int, int]]]]]] ..

static from_aligned_structures(lhs_structure, rhs_structure, barriers=None, elementary_step_id=None)[source]

Generates a reaction template from two aligned structures.

This function expects the atoms within both given structures to mach. That is that the atoms have consistent indices.

Parameters:
lhs_structurescine_utilities.AtomCollection

The left-hand side molecules in a single structure.

rhs_structurescine_utilities.AtomCollection

The right-hand side molecules in a single structure.

barriersOptional[Tuple[float, float]], optional

The lowest known barriers in kJ/mol for forward. ('lhs' -> 'rhs') and backward ('rhs' -> 'lhs') reaction.

elementary_step_idOptional[str], optional

The ID of the/an analyzed elementary step in a SCINE Database.

Returns:
Optional[ReactionTemplate]

The generated reaction template.

rtype:

Optional[ReactionTemplate] ..

static from_trajectory(trajectory, elementary_step_id=None)[source]

Generates a reaction template from a given molecular trajectory.

This function will try to extract a transition state energy and reaction barriers from the trajectory. If this behaviour is not wanted make sure that the trajectory does not contain energies.

Parameters:
trajectoryscine_utilities.MolecularTrajectory

The molecular trajectory.

elementary_step_idOptional[str], optional

The ID of the/an analyzed elementary step in a SCINE Database.

Returns:
Optional[ReactionTemplate]

The generated reaction template.

rtype:

Optional[ReactionTemplate] ..

static from_trajectory_spline(spline, reverse=False, barriers=None, elementary_step_id=None)[source]

Generates a reaction template from a molecular trajectory given as a spline interpolation.

This function will not extract energy formation from the spline.

Parameters:
splinescine_utilities.bsplines.TrajectorySpline

The molecular trajectory.

reversebool

If true will reverse the spline before analysis, False by default.

barriersOptional[Tuple[float, float]], optional

The lowest known barriers in kJ/mol for forward. ('lhs' -> 'rhs') and backward ('rhs' -> 'lhs') reaction.

elementary_step_idOptional[str], optional

The ID of the/an analyzed elementary step in a SCINE Database.

Returns:
Optional[ReactionTemplate]

The generated reaction template.

rtype:

Optional[ReactionTemplate] ..

static generate_networkx_graph_from_template(template)[source]

Generates a graph representation of the reaction template.

Parameters:
templateReactionTemplate

The template to be converted.

Returns:
Tuple[networkx.Graph, dict]

Contains the networkx.Graph and a mapping between graph nodes atoms in the template.

rtype:

Tuple[Graph, dict] ..

get_networkx_graph()[source]

Get the graph representation of this reaction template.

Returns:
Tuple[networkx.Graph, dict]

Contains the networkx.Graph and a mapping between graph nodes atoms in the template.

rtype:

Tuple[Graph, dict] ..

get_uuid()[source]
Returns:
str

A unique string identifying this template.

rtype:

str ..

strict_equal(other)[source]

At the moment this function is identical to the basic __eq__ operator, in future versions it may contain additional checks.

Parameters:
otherReactionTemplate

The other reaction template to be checked.

Returns:
bool

True if the two templates are identical.

rtype:

bool ..

update(other)[source]

If the two templates are identical (isomorphic minimal template graphs), the data of the argument will be subsumed into the base template.

This means shape information, recorded energy barriers and known elementary steps will be combined.

Parameters:
otherReactionTemplate

The other reaction template to update this with.

:rtype: :py:obj:`None`
class scine_art.reaction_template.SideContainer(lhs, rhs)[source]

A dictionary like container with the keys lhs and rhs indicating the two sides of a reaction. For properties encoding mappings from one side to the other, the given key gives the side that if mapped from.

Examples

>>> counters: SideContainer[int] = ...
>>> print(f'Number of molecules on the left-hand side: {counters["lhs"]}')
>>> print(f'Number of molecules on the right-hand side: {counters["rhs"]}')
class scine_art.reaction_template.TemplateTypeContainer(minimal=None, minimal_shell=None, fragment=None, fragment_shell=None)[source]

A dictionary like container for different template sub-types. Depending on the requested types upon analysis the keys ['minimal', 'minimal_shell', 'fragment', 'fragment_shell'], may be populated with data.

Reaction Template Database

class scine_art.database.ReactionTemplateDatabase[source]

A database of reaction templates.

add_scine_database(host='localhost', port=27017, name='default', energy_cutoff=300, model=None, gibbs_free_energy=True)[source]

Parse an entire SCINE Database adding and deduplicating all templates into this database.

Parameters:
hoststr, optional

The IP or hostname of the database server, by default ‘localhost’

portint, optional

The port of the database server, by default 27017

namestr, optional

The name of the database on the server, by default ‘default’

energy_cutofffloat, optional

Only store templates where at least one recorded barrier of the reaction is below the given threshold (in kJ/mol), by default 300 kJ/mol.

modelOptional[scine_database.Model], optional

The scine_database.Model to be allowed for energy evaluations, by default None is given, data of all models in the SCINE Database are allowed.

gibbs_free_energybool, optional

If true will only consider Gibbs free energies as valid energies for evaluation, by default True

:rtype: :py:obj:`None`
add_template(new_template)[source]

Add a new template to the database.

Will deduplicate the reaction template against stored ones.

Parameters:
new_templateReactionTemplate

The new reaction template

Returns:
str

Returns the UUID of the now stored template, if the template was a duplicate and thus subsumed, the UUID of the existing match is returned.

rtype:

str ..

append_file(path='.rtdb.pickle.obj')[source]

Loads a stored (pickled) file of templates and adds it to the existing unique reaction templates, keeps current template and adds data from the file.

Parameters:
pathstr, optional

The path to a stored list of reaction templates, by default '.rtdb.pickle.obj'

:rtype: :py:obj:`None`
dump_svg_representation(path, index=None)[source]

Generate SVG representations of all templates.

Parameters:
pathstr

The path (folder) to generate the SVG files in.

indexOptional[int], optional

If given only generate an SVG for the template with the given index, by default None

:rtype: :py:obj:`None`
find_matching_templates(molecules, energy_cutoff=300, enforce_atom_shapes=True)[source]

Finds all matching reaction templates for a given set of molecules.

All given molecules must be used in the matched reaction template exactly once. None may be unused or used twice.

A single template can match multiple atom combinations in the given atoms, hence multiple matches can be returned per reaction template.

Parameters:
moleculesList[masm.Molecule]

The molecules to apply the template to.

energy_cutofffloat, optional

Only apply the template if a barrier in the trialed direction of less than this value has been reported, by default 300 (kJ/mol)

enforce_atom_shapesbool, optional

If true only allow atoms with the same coordination sphere shapes to be considered matching, by default True.

Returns:
Optional[List[Match]]

If any matches exist, returns a list of all possible matches.

rtype:

Optional[List[Dict[str, List[Tuple[Tuple[int, int], Tuple[int, int]]]]]] ..

get_template(uuid)[source]

Getter for a single reaction template.

Parameters:
uuidstr

The unique ID of a reaction template.

Returns:
Optional[ReactionTemplate]

The reaction template with the given unique ID if present.

rtype:

Optional[ReactionTemplate] ..

iterate_templates()[source]

Iterator through all stored reaction templates.

Yields:
Generator[ReactionTemplate, None, None]

Iterator of all reaction templates in the database.

rtype:

Generator[ReactionTemplate, None, None] ..

load(path='.rtdb.pickle.obj')[source]

Loads a stored (pickled) file of templates and replaces existing data.

Parameters:
pathstr, optional

The path to a stored list of reaction templates, by default '.rtdb.pickle.obj'

:rtype: :py:obj:`None`
save(path='.rtdb.pickle.obj')[source]

Saves the currently stored unique reaction templates into a file.

Parameters:
pathstr, optional

The path to store a list of reaction templates, by default '.rtdb.pickle.obj'

:rtype: :py:obj:`None`
template_count()[source]

Getter for the number of reaction templates int he database.

Returns:
int

Number of unique reaction templates in the database.

rtype:

int ..

Molecule Comparisons

scine_art.molecules.get_atom_mapping(m1, m2)[source]

Maps atoms between two identical molecules.

Parameters:
m1scine_molassembler.Molecule

The first molecule.

m2scine_molassembler.Molecule

The second molecule.

Returns:
List[int]

Mapping of atom indices of the first molecule to those of the second molecule.

Raises:
RuntimeError

If molecules are not identical.

Examples

>>> mapping = get_atom_mapping(m1, m2)
>>> for i in range(m1.graph.V):
>>>     assert m1.graph.element_type() == m2.graph.element_type(mapping[i])
Return type:

List[int]

scine_art.molecules.intersection(m1, m2)[source]

Identifies all intersecting fragments in two molecules.

Parameters:
m1scine_molassembler.Molecule

The first molecule.

m2scine_molassembler.Molecule

The second molecule.

Returns:
Tuple[List[scine_molassembler.Molecule], List[Tuple[List[List[int]], List[List[int]]]]]

A list of unique molecular fragments that are present in both molecules. Additionally, for each fragment all mappings of atom indices to each of the two molecules are returned.

Examples

>>> for fragment, mappings in intersection(m1, m2):
>>>     print(mappings[0]) #  all mappings of fragment in m1
>>>     print(mappings[1]) #  all mappings of fragment in m2
Return type:

Tuple[List[Molecule], List[Tuple[List[List[int]], List[List[int]]]]]

scine_art.molecules.is_single_subgraph(ids, graph)[source]

Checks if the set of indices given describes a single, connected subgraph/fragment of the given molecular graph.

Parameters:
idsList[int]

List of atom indices.

graphscine_molassembler.Graph

Molecular graph.

Returns:
bool

True if the indices describe one connected subgraph. False otherwise.

rtype:

bool ..

scine_art.molecules.matching_atom_ids(m1, m2)[source]

Generates sets of atom indices of all matching molecular fragments.

Parameters:
m1scine_molassembler.Molecule

The first molecule.

m2scine_molassembler.Molecule

The second molecule.

Returns:
List[Tuple[Set[int], Set[int]]]

A list of pairs of atom index sets. The sets represent atom indices in the respective molecules. Each pair represents one matching fragment. Fragments are non-unique, as a single fragment of one molecule may match multiple in the other.

rtype:

List[Tuple[Set[int], Set[int]]] ..

scine_art.molecules.mol_from_subgraph_indices(subgraph_indices, molecule)[source]

Cuts down molecule to fragment identified by atom indices.

Additionally produces direct atom mapping between the initial molecule and the resulting fragment.

Parameters:
subgraph_indicesList[int]

Atom indices of a connected subgraph in the given molecule.

moleculescine_molassembler.Molecule

The molecule to reduce into a fragment.

Returns:
Tuple[Optional[scine_molassembler.Molecule], List[int]]

The molecular fragment and its atom mapping towards the original molecule. The molecule returned will be None if the initial index list was empty.

rtype:

Tuple[Optional[Molecule], List[int]] ..

scine_art.molecules.set_difference(m1, m2)[source]

_summary_

Parameters:
m1scine_molassembler.Molecule

_description_

m2scine_molassembler.Molecule

_description_

Returns:
Tuple[List[scine_molassembler.Molecule], List[List[List[int]]]]

_description_

rtype:

Tuple[List[Molecule], List[List[List[int]]]] ..

scine_art.molecules.sort_indices_by_subgraph(ids, graph)[source]

Groups the given indices by subgraphs.

All indices that are directly connected to one another (connected in a single subgraph) are sorted together. All indices will be returned, in a limiting case each index will be part of its own subgraph.

Parameters:
idsList[int]

Indices to be sorted.

graphscine_molassembler.Graph

Molecular graph to sort the indices by.

Returns:
List[List[int]]

A list of subgraphs, each given as list of indices.

rtype:

List[List[int]] ..

scine_art.molecules.symmetric_difference(m1, m2)[source]

_summary_

Parameters:
m1scine_molassembler.Molecule

_description_

m2scine_molassembler.Molecule

_description_

Returns:
Tuple[List[scine_molassembler.Molecule], List[Tuple[List[List[int]], List[List[int]]]]]

_description_

rtype:

Tuple[List[Molecule], List[Tuple[List[List[int]], List[List[int]]]]] ..

IO

scine_art.io.load_file(path)[source]

Load an atom coordinate file as scine_molassembler.Molecule s.

If the fle format does not provide bond information, then the bonds will be inferred using distance criteria based on covalent radii. For more information on there bond inference see the scine_utilities.BondDetector.

Parameters:
pathstr

Path to the structure file.

Returns:
List[scine_molassembler.Molecule]

The resulting interpreted molecules.

scine_art.io.load_spline_from_trajectory(path)[source]

Read a trajectory from an .xyz file.

Expects energies to be given as plain floating point numbers in the title field of each structure

Parameters:
pathstr

Path to the trajectory file (.xyz)

Returns:
utils.bsplines.TrajectorySpline

The spline interpolation of the read trajectory.

rtype:

TrajectorySpline ..

scine_art.io.write_molecule_to_svg(path, mol)[source]

Dump a single molecular graph to disk in SVG format.

Parameters:
pathstr

Path to the SVG file to be generated.

molscine_molassembler.Molecule

The molecule to be dumped to disk as an SVG file.

:rtype: :py:obj:`None`

Experimental: Atom Matching

scine_art.experimental.map_reaction_from_molecules_cached(lhs, rhs, known_atom_mappings=None)[source]

Tries to match each atom in the molecules of the left-hand side to exactly one atom on the right-hand side.

Atom type counts must match.

Parameters:
lhsList[masm.Molecule]

A list of molecules on the left-hand side of the matching.

rhsList[masm.Molecule]

A list of molecules on the right-hand side of the matching.

known_atom_mappingsOptional[List[Tuple[int, int]]]

A list of atom mappings that have to be respected. The atoms are indexed on a continuous scale in order of the atoms in the molecules given.

Returns:
Tuple[List[List[Tuple[int, int]]], List[List[Tuple[int, int]]]]

A tuple atom matches. The first entry matching from left- to right-hand side the second from right- to left-hand side.

Raises:
RuntimeError

If algorithm fails or miss matches appear.

rtype:

Tuple[List[List[Tuple[int, int]]], List[List[Tuple[int, int]]]] ..

scine_art.experimental.map_reaction_from_molecules_direct(lhs, rhs, known_atom_mappings=None)[source]

Tries to match each atom in the molecules of the left-hand side to exactly one atom on the right-hand side.

Atom type counts must match.

Parameters:
lhsList[masm.Molecule]

A list of molecules on the left-hand side of the matching.

rhsList[masm.Molecule]

A list of molecules on the right-hand side of the matching.

known_atom_mappingsOptional[List[Tuple[int, int]]]

A list of atom mappings that have to be respected. The atoms are indexed on a continuous scale in order of the atoms in the molecules given.

Returns:
Tuple[List[List[Tuple[int, int]]], List[List[Tuple[int, int]]]]

A tuple atom matches. The first entry matching from left- to right-hand side the second from right- to left-hand side.

Raises:
RuntimeError

If algorithm fails or miss matches appear.

rtype:

Tuple[List[List[Tuple[int, int]]], List[List[Tuple[int, int]]]] ..