rmgpy.molecule.Molecule¶
- class rmgpy.molecule.Molecule¶
A representation of a molecular structure using a graph data type, extending the
Graph
class. Attributes are:Attribute
Type
Description
atoms
list
A list of Atom objects in the molecule
symmetry_number
float
The (estimated) external + internal symmetry number of the molecule, modified for chirality
multiplicity
int
The multiplicity of this species, multiplicity = 2*total_spin+1
reactive
bool
True
(by default) if the molecule participates in reaction families.It is set to
False
by the filtration functions if a non representative resonance structure was generated by a template reaction
props
dict
A list of properties describing the state of the molecule.
inchi
str
A string representation of the molecule in InChI
smiles
str
A string representation of the molecule in SMILES
fingerprint
str
A representation for fast comparison, set as molecular formula
metal
str
The metal of the metal surface the molecule is associated with
facet
str
The facet of the metal surface the molecule is associated with
A new molecule object can be easily instantiated by passing the smiles or inchi string representing the molecular structure.
- add_atom(atom)¶
Add an atom to the graph. The atom is initialized with no bonds.
- add_bond(bond)¶
Add a bond to the graph as an edge connecting the two atoms atom1 and atom2.
- add_edge(edge)¶
Add an edge to the graph. The two vertices in the edge must already exist in the graph, or a
ValueError
is raised.
- add_vertex(vertex)¶
Add a vertex to the graph. The vertex is initialized with no edges.
- assign_atom_ids()¶
Assigns an index to every atom in the molecule for tracking purposes. Uses entire range of cython’s integer values to reduce chance of duplicates
- atom_ids_valid()¶
Checks to see if the atom IDs are valid in this structure
- atoms¶
List of atoms contained in the current molecule.
Renames the inherited vertices attribute of
Graph
.
- calculate_cp0()¶
Return the value of the heat capacity at zero temperature in J/mol*K.
- calculate_cpinf()¶
Return the value of the heat capacity at infinite temperature in J/mol*K.
- calculate_symmetry_number()¶
Return the symmetry number for the structure. The symmetry number includes both external and internal modes.
- clear_labeled_atoms()¶
Remove the labels from all atoms in the molecule.
- connect_the_dots(critical_distance_factor, raise_atomtype_exception)¶
Delete all bonds, and set them again based on the Atoms’ coords. Does not detect bond type.
- contains_labeled_atom(label)¶
Return
True
if the molecule contains an atom with the label label andFalse
otherwise.
- contains_surface_site()¶
Returns
True
iff the molecule contains an ‘X’ surface site.
- copy(deep)¶
Create a copy of the current graph. If deep is
True
, a deep copy is made: copies of the vertices and edges are used in the new graph. If deep isFalse
or not specified, a shallow copy is made: the original vertices and edges are used in the new graph.
- copy_and_map()¶
Create a deep copy of the current graph, and return the dict ‘mapping’. Method was modified from Graph.copy() method
- count_aromatic_rings()¶
Count the number of aromatic rings in the current molecule, as determined by the benzene bond type. This is purely dependent on representation and is unrelated to the actual aromaticity of the molecule.
Returns an integer corresponding to the number or aromatic rings.
- count_internal_rotors()¶
Determine the number of internal rotors in the structure. Any single bond not in a cycle and between two atoms that also have other bonds is considered to be a pivot of an internal rotor.
- delete_hydrogens()¶
Irreversibly delete all non-labeled hydrogens without updating connectivity values. If there’s nothing but hydrogens, it does nothing. It destroys information; be careful with it.
- draw(path)¶
Generate a pictorial representation of the chemical graph using the
draw
module. Use path to specify the file to save the generated image to; the image type is automatically determined by extension. Valid extensions are.png
,.svg
,.pdf
, and.ps
; of these, the first is a raster format and the remainder are vector formats.
- enumerate_bonds()¶
Count the number of each type of bond (e.g. ‘C-H’, ‘C=C’) present in the molecule :return: dictionary, with bond strings as keys and counts as values
- find_h_bonds()¶
generates a list of (new-existing H bonds ignored) possible Hbond coordinates [(i1,j1),(i2,j2),…] where i and j values correspond to the indexes of the atoms involved, Hbonds are allowed if they meet the following constraints:
between a H and [O,N] atoms
the hydrogen is covalently bonded to an O or N
the Hydrogen bond must complete a ring with at least 5 members
An atom can only be hydrogen bonded to one other atom
- find_isomorphism(other, initial_map, save_order, strict)¶
Returns
True
if other is isomorphic andFalse
otherwise, and the matching mapping. The initialMap attribute can be used to specify a required mapping from self to other (i.e. the atoms of self are the keys, while the atoms of other are the values). The returned mapping also uses the atoms of self for the keys and the atoms of other for the values. The other parameter must be aMolecule
object, or aTypeError
is raised.- Parameters:
initial_map (dict, optional) – initial atom mapping to use
save_order (bool, optional) – if
True
, reset atom order after performing atom isomorphismstrict (bool, optional) – if
False
, perform isomorphism ignoring electrons
- find_subgraph_isomorphisms(other, initial_map, save_order)¶
Returns
True
if other is subgraph isomorphic andFalse
otherwise. Also returns the lists all of valid mappings. The initial_map attribute can be used to specify a required mapping from self to other (i.e. the atoms of self are the keys, while the atoms of other are the values). The returned mappings also use the atoms of self for the keys and the atoms of other for the values. The other parameter must be aGroup
object, or aTypeError
is raised.
- fingerprint¶
Fingerprint used to accelerate graph isomorphism comparisons with other molecules. The fingerprint is a short string containing a summary of selected information about the molecule. Two fingerprint strings matching is a necessary (but not sufficient) condition for the associated molecules to be isomorphic.
Use an expanded molecular formula to also enable sorting.
- from_adjacency_list(adjlist, saturate_h, raise_atomtype_exception, raise_charge_exception, check_consistency)¶
Convert a string adjacency list adjlist to a molecular structure. Skips the first line (assuming it’s a label) unless withLabel is
False
.
- from_augmented_inchi(aug_inchi, raise_atomtype_exception)¶
Convert an Augmented InChI string aug_inchi to a molecular structure.
- from_inchi(inchistr, backend, raise_atomtype_exception)¶
Convert an InChI string inchistr to a molecular structure.
RDKit and Open Babel are the two backends used in RMG. It is possible to use a single backend or try different backends in sequence. The available options for the
backend
argument: ‘openbabel-first’(default), ‘rdkit-first’, ‘rdkit’, or ‘openbabel’.
- from_smarts(smartsstr, raise_atomtype_exception)¶
Convert a SMARTS string smartsstr to a molecular structure. Uses RDKit to perform the conversion. This Kekulizes everything, removing all aromatic atom types.
- from_smiles(smilesstr, backend, raise_atomtype_exception)¶
Convert a SMILES string smilesstr to a molecular structure.
RDKit and Open Babel are the two backends used in RMG. It is possible to use a single backend or try different backends in sequence. The available options for the
backend
argument: ‘openbabel-first’(default), ‘rdkit-first’, ‘rdkit’, or ‘openbabel’.
- from_xyz(atomic_nums, coordinates, critical_distance_factor, raise_atomtype_exception)¶
Create an RMG molecule from a list of coordinates and a corresponding list of atomic numbers. These are typically received from CCLib and the molecule is sent to ConnectTheDots so will only contain single bonds.
- generate_h_bonded_structures()¶
generates a list of Hbonded molecular structures in addition to the constraints on Hydrogen bonds applied in the find_H_Bonds function the generated structures are constrained to:
An atom can only be hydrogen bonded to one other atom
Only two H-bonds can exist in a given molecule
the second is done to avoid explosive growth in the number of structures as without this constraint the number of possible structures grows 2^n where n is the number of possible H-bonds
- generate_resonance_structures(keep_isomorphic, filter_structures, save_order)¶
Returns a list of resonance structures of the molecule.
- get_adatoms()¶
Get a list of adatoms in the molecule. :returns: A list containing the adatoms in the molecule :rtype: List(Atom)
- get_all_cycles(starting_vertex)¶
Given a starting vertex, returns a list of all the cycles containing that vertex.
This function returns a duplicate of each cycle because [0,1,2,3] is counted as separate from [0,3,2,1]
- get_all_cycles_of_size(size)¶
Return a list of the all non-duplicate rings with length ‘size’. The algorithm implements was adapted from a description by Fan, Panaye, Doucet, and Barbu (doi: 10.1021/ci00015a002)
B. T. Fan, A. Panaye, J. P. Doucet, and A. Barbu. “Ring Perception: A New Algorithm for Directly Finding the Smallest Set of Smallest Rings from a Connection Table.” J. Chem. Inf. Comput. Sci. 33, p. 657-662 (1993).
- get_all_cyclic_vertices()¶
Returns all vertices belonging to one or more cycles.
- get_all_edges()¶
Returns a list of all edges in the graph.
- get_all_labeled_atoms()¶
Return the labeled atoms as a
dict
with the keys being the labels and the values the atoms themselves. If two or more atoms have the same label, the value is converted to a list of these atoms.
- get_all_polycyclic_vertices()¶
Return all vertices belonging to two or more cycles, fused or spirocyclic.
- get_all_simple_cycles_of_size(size)¶
Return a list of all non-duplicate monocyclic rings with length ‘size’.
Naive approach by eliminating polycyclic rings that are returned by
getAllCyclicsOfSize
.
- get_aromatic_rings(rings, save_order)¶
Returns all aromatic rings as a list of atoms and a list of bonds.
Identifies rings using Graph.get_smallest_set_of_smallest_rings(), then uses RDKit to perceive aromaticity. RDKit uses an atom-based pi-electron counting algorithm to check aromaticity based on Huckel’s Rule. Therefore, this method identifies “true” aromaticity, rather than simply the RMG bond type.
The method currently restricts aromaticity to six-membered carbon-only rings. This is a limitation imposed by RMG, and not by RDKit.
By default, the atom order will be sorted to get consistent results from different runs. The atom order can be saved when dealing with problems that are sensitive to the atom map.
- get_bond(atom1, atom2)¶
Returns the bond connecting atoms atom1 and atom2.
- get_bonds(atom)¶
Return a dictionary of the bonds involving the specified atom.
- get_charge_span()¶
Iterate through the atoms in the structure and calculate the charge span on the overall molecule. The charge span is a measure of the number of charge separations in a molecule.
- get_desorbed_molecules()¶
Get a list of desorbed molecules by desorbing the molecule from the surface.
Returns a list of Molecules. Each molecule’s atoms will be labeled corresponding to the bond order with the surface:
*1
- single bond*2
- double bond*3
- triple bond*4
- quadruple bond
- get_deterministic_sssr()¶
Modified Graph method get_smallest_set_of_smallest_rings by sorting calculated cycles by short length and then high atomic number instead of just short length (for cases where multiple cycles with same length are found, get_smallest_set_of_smallest_rings outputs non-determinstically).
For instance, molecule with this smiles: C1CC2C3CSC(CO3)C2C1, will have non-deterministic output from get_smallest_set_of_smallest_rings, which leads to non-deterministic bicyclic decomposition. Using this new method can effectively prevent this situation.
Important Note: This method returns an incorrect set of SSSR in certain molecules (such as cubane). It is recommended to use the main Graph.get_smallest_set_of_smallest_rings method in new applications. Alternatively, consider using Graph.get_relevant_cycles for deterministic output.
In future development, this method should ideally be replaced by some method to select a deterministic set of SSSR from the set of Relevant Cycles, as that would be a more robust solution.
- get_disparate_cycles()¶
Get all disjoint monocyclic and polycyclic cycle clusters in the molecule. Takes the RC and recursively merges all cycles which share vertices.
Returns: monocyclic_cycles, polycyclic_cycles
- get_edge(vertex1, vertex2)¶
Returns the edge connecting vertices vertex1 and vertex2.
- get_edges(vertex)¶
Return a dictionary of the edges involving the specified vertex.
- get_edges_in_cycle(vertices, sort)¶
For a given list of atoms comprising a ring, return the set of bonds connecting them, in order around the ring.
If sort=True, then sort the vertices to match their connectivity. Otherwise, assumes that they are already sorted, which is true for cycles returned by get_relevant_cycles or get_smallest_set_of_smallest_rings.
- get_element_count()¶
Returns the element count for the molecule as a dictionary.
- get_formula()¶
Return the molecular formula for the molecule.
- get_labeled_atoms(label)¶
Return the atoms in the molecule that are labeled.
- get_largest_ring(vertex)¶
returns the largest ring containing vertex. This is typically useful for finding the longest path in a polycyclic ring, since the polycyclic rings returned from get_polycycles are not necessarily in order in the ring structure.
- get_max_cycle_overlap()¶
Return the maximum number of vertices that are shared between any two cycles in the graph. For example, if there are only disparate monocycles or no cycles, the maximum overlap is zero; if there are “spiro” cycles, it is one; if there are “fused” cycles, it is two; and if there are “bridged” cycles, it is three.
- get_molecular_weight()¶
Return the molecular weight of the molecule in kg/mol.
- get_monocycles()¶
Return a list of cycles that are monocyclic.
- get_net_charge()¶
Iterate through the atoms in the structure and calculate the net charge on the overall molecule.
- get_nth_neighbor(starting_atoms, distance_list, ignore_list, n)¶
Recursively get the Nth nonHydrogen neighbors of the starting_atoms, and return them in a list. starting_atoms is a list of :class:Atom for which we will get the nth neighbor. distance_list is a list of integers, corresponding to the desired neighbor distances. ignore_list is a list of :class:Atom that have been counted in (n-1)th neighbor, and will not be returned. n is an integer, corresponding to the distance to be calculated in the current iteration.
- get_num_atoms(element)¶
Return the number of atoms in molecule. If element is given, ie. “H” or “C”, the number of atoms of that element is returned.
- get_polycycles()¶
Return a list of cycles that are polycyclic. In other words, merge the cycles which are fused or spirocyclic into a single polycyclic cycle, and return only those cycles. Cycles which are not polycyclic are not returned.
- get_radical_atoms()¶
Return the atoms in the molecule that have unpaired electrons.
- get_radical_count()¶
Return the total number of radical electrons on all atoms in the molecule. In this function, monoradical atoms count as one, biradicals count as two, etc.
- get_relevant_cycles()¶
Returns the set of relevant cycles as a list of lists. Uses RingDecomposerLib for ring perception.
Kolodzik, A.; Urbaczek, S.; Rarey, M. Unique Ring Families: A Chemically Meaningful Description of Molecular Ring Topologies. J. Chem. Inf. Model., 2012, 52 (8), pp 2013-2021
Flachsenberg, F.; Andresen, N.; Rarey, M. RingDecomposerLib: An Open-Source Implementation of Unique Ring Families and Other Cycle Bases. J. Chem. Inf. Model., 2017, 57 (2), pp 122-126
- get_singlet_carbene_count()¶
Return the total number of singlet carbenes (lone pair on a carbon atom) in the molecule. Counts the number of carbon atoms with a lone pair. In the case of [C] with two lone pairs, this method will return 1.
- get_smallest_set_of_smallest_rings()¶
Returns the smallest set of smallest rings as a list of lists. Uses RingDecomposerLib for ring perception.
Kolodzik, A.; Urbaczek, S.; Rarey, M. Unique Ring Families: A Chemically Meaningful Description of Molecular Ring Topologies. J. Chem. Inf. Model., 2012, 52 (8), pp 2013-2021
Flachsenberg, F.; Andresen, N.; Rarey, M. RingDecomposerLib: An Open-Source Implementation of Unique Ring Families and Other Cycle Bases. J. Chem. Inf. Model., 2017, 57 (2), pp 122-126
- get_surface_sites()¶
Get a list of surface site atoms in the molecule. :returns: A list containing the surface site atoms in the molecule :rtype: List(Atom)
- get_symmetry_number()¶
Returns the symmetry number of Molecule. First checks whether the value is stored as an attribute of Molecule. If not, it calls the calculate_symmetry_number method.
- get_url()¶
Get a URL to the molecule’s info page on the RMG website.
- has_atom(atom)¶
Returns
True
if atom is an atom in the graph, orFalse
if not.
- has_bond(atom1, atom2)¶
Returns
True
if atoms atom1 and atom2 are connected by an bond, orFalse
if not.
- has_edge(vertex1, vertex2)¶
Returns
True
if vertices vertex1 and vertex2 are connected by an edge, orFalse
if not.
- has_halogen()¶
Return
True
if the molecule contains at least one halogen (F, Cl, Br, or I), orFalse
otherwise.
- has_lone_pairs()¶
Return
True
if the molecule contains at least one lone electron pair, orFalse
otherwise.
- has_vertex(vertex)¶
Returns
True
if vertex is a vertex in the graph, orFalse
if not.
- identify_ring_membership()¶
Performs ring perception and saves ring membership information to the Atom.props attribute.
- inchi¶
InChI string for this molecule. Read-only.
- is_aromatic()¶
Returns
True
if the molecule is aromatic, orFalse
if not. Iterates over the SSSR’s and searches for rings that consist solely of Cb atoms. Assumes that aromatic rings always consist of 6 atoms. In cases of naphthalene, where a 6 + 4 aromatic system exists, there will be at least one 6 membered aromatic ring so this algorithm will not fail for fused aromatic rings.
- is_aryl_radical(aromatic_rings, save_order)¶
Return
True
if the molecule only contains aryl radicals, i.e., radical on an aromatic ring, orFalse
otherwise. If noaromatic_rings
provided, aromatic rings will be searched in-place, and this process may involve atom order change by default. Setsave_order
toTrue
to force the atom order unchanged.
- is_atom_in_cycle(atom)¶
Return :data:
True
ifatom
is in one or more cycles in the structure, and :data:False
if not.
- is_bond_in_cycle(bond)¶
Return :data:
True
if the bond between atomsatom1
andatom2
is in one or more cycles in the graph, or :data:False
if not.
- is_cyclic()¶
Return
True
if one or more cycles are present in the graph orFalse
otherwise.
- is_edge_in_cycle(edge)¶
Return
True
if the edge between vertices vertex1 and vertex2 is in one or more cycles in the graph, orFalse
if not.
- is_heterocyclic()¶
Returns
True
if the molecule is heterocyclic, orFalse
if not.
- is_identical(other, strict)¶
Performs isomorphism checking, with the added constraint that atom IDs must match.
Primary use case is tracking atoms in reactions for reaction degeneracy determination.
Returns
True
if two graphs are identical andFalse
otherwise.If
strict=False
, performs the check ignoring electrons and resonance structures.
- is_isomorphic(other, initial_map, generate_initial_map, save_order, strict)¶
Returns
True
if two graphs are isomorphic andFalse
otherwise. The initialMap attribute can be used to specify a required mapping from self to other (i.e. the atoms of self are the keys, while the atoms of other are the values). The other parameter must be aMolecule
object, or aTypeError
is raised. Also ensures multiplicities are also equal.- Parameters:
initial_map (dict, optional) – initial atom mapping to use
generate_initial_map (bool, optional) – if
True
, initialize map by pairing atoms with same labelssave_order (bool, optional) – if
True
, reset atom order after performing atom isomorphismstrict (bool, optional) – if
False
, perform isomorphism ignoring electrons
- is_linear()¶
Return
True
if the structure is linear andFalse
otherwise.
- is_mapping_valid(other, mapping, equivalent, strict)¶
Check that a proposed mapping of vertices from self to other is valid by checking that the vertices and edges involved in the mapping are mutually equivalent. If equivalent is
True
it checks if atoms and edges are equivalent, ifFalse
it checks if they are specific cases of each other. If strict isTrue
, electrons and bond orders are considered, and ignored ifFalse
.
- is_radical()¶
Return
True
if the molecule contains at least one radical electron, orFalse
otherwise.
- is_subgraph_isomorphic(other, initial_map, generate_initial_map, save_order)¶
Returns
True
if other is subgraph isomorphic andFalse
otherwise. The initial_map attribute can be used to specify a required mapping from self to other (i.e. the atoms of self are the keys, while the atoms of other are the values). The other parameter must be aGroup
object, or aTypeError
is raised.
- is_surface_site()¶
Returns
True
iff the molecule is nothing but a surface site ‘X’.
- is_vertex_in_cycle(vertex)¶
Return
True
if the given vertex is contained in one or more cycles in the graph, orFalse
if not.
- kekulize()¶
Kekulizes an aromatic molecule.
- merge(other)¶
Merge two molecules so as to store them in a single
Molecule
object. The mergedMolecule
object is returned.
- number_of_surface_sites()¶
Returns the number of surface sites in the molecule. e.g. 2 for a bidentate adsorbate
- remove_atom(atom)¶
Remove atom and all bonds associated with it from the graph. Does not remove atoms that no longer have any bonds as a result of this removal.
- remove_bond(bond)¶
Remove the bond between atoms atom1 and atom2 from the graph. Does not remove atoms that no longer have any bonds as a result of this removal.
- remove_edge(edge)¶
Remove the specified edge from the graph. Does not remove vertices that no longer have any edges as a result of this removal.
- remove_h_bonds()¶
removes any present hydrogen bonds from the molecule
- remove_van_der_waals_bonds()¶
Remove all van der Waals bonds.
- remove_vertex(vertex)¶
Remove vertex and all edges associated with it from the graph. Does not remove vertices that no longer have any edges as a result of this removal.
- replace_halogen_with_hydrogen(raise_atomtype_exception)¶
Replace all halogens in a molecule with hydrogen atoms. Changes self molecule object.
- reset_connectivity_values()¶
Reset any cached connectivity information. Call this method when you have modified the graph.
- restore_vertex_order()¶
reorder the vertices to what they were before sorting if you saved the order
- saturate_radicals(raise_atomtype_exception)¶
Saturate the molecule by replacing all radicals with bonds to hydrogen atoms. Changes self molecule object.
- saturate_unfilled_valence(update)¶
Saturate the molecule by adding H atoms to any unfilled valence
- smiles¶
SMILES string for this molecule. Read-only.
- sort_atoms()¶
Sort the atoms in the graph. This can make certain operations, e.g. the isomorphism functions, much more efficient.
This function orders atoms using several attributes in atom.getDescriptor(). Currently it sorts by placing heaviest atoms first and hydrogen atoms last. Placing hydrogens last during sorting ensures that functions with hydrogen removal work properly.
- sort_cyclic_vertices(vertices)¶
Given a list of vertices comprising a cycle, sort them such that adjacent entries in the list are connected to each other. Warning: Assumes that the cycle is elementary, ie. no bridges.
- sort_vertices(save_order)¶
Sort the vertices in the graph. This can make certain operations, e.g. the isomorphism functions, much more efficient.
- sorting_key¶
Returns a sorting key for comparing Molecule objects. Read-only
- split()¶
Convert a single
Molecule
object containing two or more unconnected molecules into separate class:Molecule objects.
- to_adjacency_list(label, remove_h, remove_lone_pairs, old_style)¶
Convert the molecular structure to a string adjacency list.
- to_augmented_inchi(backend)¶
Adds an extra layer to the InChI denoting the multiplicity of the molecule.
Separate layer with a forward slash character.
RDKit and Open Babel are the two backends used in RMG. It is possible to use a single backend or try different backends in sequence. The available options for the
backend
argument: ‘rdkit-first’(default), ‘openbabel-first’, ‘rdkit’, or ‘openbabel’.
- to_augmented_inchi_key(backend)¶
Adds an extra layer to the InChIKey denoting the multiplicity of the molecule.
Simply append the multiplicity string, do not separate by a character like forward slash.
RDKit and Open Babel are the two backends used in RMG. It is possible to use a single backend or try different backends in sequence. The available options for the
backend
argument: ‘rdkit-first’(default), ‘openbabel-first’, ‘rdkit’, or ‘openbabel’.
- to_group()¶
This method converts a list of atoms in a Molecule to a Group object.
- to_inchi(backend)¶
Convert a molecular structure to an InChI string. Uses RDKit to perform the conversion. Perceives aromaticity.
or
Convert a molecular structure to an InChI string. Uses OpenBabel to perform the conversion.
It is possible to use a single backend or try different backends in sequence. The available options for the
backend
argument: ‘rdkit-first’(default), ‘openbabel-first’, ‘rdkit’, or ‘openbabel’.
- to_inchi_key(backend)¶
Convert a molecular structure to an InChI Key string. Uses OpenBabel to perform the conversion.
or
Convert a molecular structure to an InChI Key string. Uses RDKit to perform the conversion.
It is possible to use a single backend or try different backends in sequence. The available options for the
backend
argument: ‘rdkit-first’(default), ‘openbabel-first’, ‘rdkit’, or ‘openbabel’.
- to_rdkit_mol(*args, **kwargs)¶
Convert a molecular structure to a RDKit rdmol object.
- to_single_bonds(raise_atomtype_exception)¶
Returns a copy of the current molecule, consisting of only single bonds.
This is useful for isomorphism comparison against something that was made via from_xyz, which does not attempt to perceive bond orders
- to_smarts()¶
Convert a molecular structure to an SMARTS string. Uses RDKit to perform the conversion. Perceives aromaticity and removes Hydrogen atoms.
- to_smiles()¶
Convert a molecular structure to an SMILES string.
If there is a Nitrogen atom present it uses OpenBabel to perform the conversion, and the SMILES may or may not be canonical.
Otherwise, it uses RDKit to perform the conversion, so it will be canonical SMILES. While converting to an RDMolecule it will perceive aromaticity and removes Hydrogen atoms.
- update(log_species, raise_atomtype_exception, sort_atoms)¶
Update the charge and atom types of atoms. Update multiplicity, and sort atoms (if
sort_atoms
isTrue
) Does not necessarily update the connectivity values (which are used in isomorphism checks) If you need that, call update_connectivity_values()
- update_atomtypes(log_species, raise_exception)¶
Iterate through the atoms in the structure, checking their atom types to ensure they are correct (i.e. accurately describe their local bond environment) and complete (i.e. are as detailed as possible).
If raise_exception is False, then the generic atomtype ‘R’ will be prescribed to any atom when get_atomtype fails. Currently used for resonance hybrid atom types.
- update_connectivity_values()¶
Update the connectivity values for each vertex in the graph. These are used to accelerate the isomorphism checking.
- update_lone_pairs()¶
Iterate through the atoms in the structure and calculate the number of lone electron pairs, assuming a neutral molecule.
- update_multiplicity()¶
Update the multiplicity of a newly formed molecule.