Input / Output¶
IO Submodule
-
class
scine_molassembler.io.
LineNotation
¶ Generates
Molecule
instances from line notations of molecules via OpenBabel, if found in the runtime path.-
__init__
(*args, **kwargs)¶ Initialize self. See help(type(self)) for accurate signature.
-
enabled
= True¶
-
static
from_canonical_smiles
(canonical_smiles: str) → scine_molassembler.Molecule¶ Construct a single
Molecule
from a canonical SMILES string
-
static
from_inchi
(inchi: str) → scine_molassembler.Molecule¶ Construct a single
Molecule
from an InChI string
-
static
from_isomeric_smiles
(isomeric_smiles: str) → scine_molassembler.Molecule¶ Construct a single
Molecule
from an isomeric SMILES string
-
-
scine_molassembler.io.
read
(filename: str) → scine_molassembler.Molecule¶ Reads a single
Molecule
from a file. Interprets the file format from its extension. Supported formats: - mol: MOLFile V2000 - xyz: XYZ file - cbor/bson/json: Serialization formats of molecules- Parameters
filename – File to read.
-
scine_molassembler.io.
split
(filename: str) → List[scine_molassembler.Molecule]¶ Reads multiple molecules from a file. Interprets the file format from its extension just like read(). Note that serializations of molecules contain only a single
Molecule
. Use read() instead.- Parameters
filename – File to read.
-
scine_molassembler.io.
write
(*args, **kwargs)¶ Overloaded function.
write(filename: str, molecule: scine_molassembler.Molecule, positions: numpy.ndarray[numpy.float64[m, 3]]) -> None
Write a
Molecule
and its positions to a file- param filename
File to write to. File format is interpreted from this parameter’s file extension.
- param molecule
Molecule
to write to file- param positions
Positions of molecule’s atoms in bohr
write(filename: str, molecule: scine_molassembler.Molecule) -> None
Write a
Molecule
serialization with the endings json/cbor/bson or a graph representation with ending dot/svg to a file.- param filename
File to write to. File format is interpreted from this parameter’s file extension
- param molecule
Molecule
to write to file
Experimental¶
- note
Functions in this module are unstable and should be used with caution. Check your results. Upon stabilization, functions will be deprecated and move to a different module.
-
scine_molassembler.io.experimental.
emit_smiles
(molecule: scine_molassembler.Molecule) → str¶ Generate a smiles string for a molecule
- Parameters
molecule – Molecule to generate smiles string for
- Returns
A (partially) normalized openSMILES-standard compliant smiles string.
- Warning
This is a lossy serialization format! The openSMILES standard does not contain stereodescriptors for shapes other than the tetrahedron, square, trigonal bipyramid and octahedron. Generated smiles containing stereocenters with other shapes will not contain stereodescriptors for these centers.
- Note
Missing normalization: Aromaticity detection in kekulized graph to aromatic atom types.
>>> biphenyl = io.experimental.from_smiles("c1ccccc1-c2ccccc2") >>> io.experimental.emit_smiles(biphenyl) 'c1ccccc1-c2ccccc2'
-
scine_molassembler.io.experimental.
from_smiles
(smiles_str: str) → scine_molassembler.Molecule¶ Parse a smiles string containing only a single molecule
The smiles parser is implemented according to the OpenSMILES spec. It supports the following features:
Isotope markers
Valence filling of the organic subset
Set local shapes from VSEPR
Ring closures (and concatenation between dot-separated components)
Stereo markers - Double bond - Tetrahedral - Square planar - Trigonal bipyramidal - Octahedral
- Parameters
smiles_str – A smiles string containing a single molecule
>>> import scine_utilities as utils >>> methane = io.experimental.from_smiles("C") >>> methane.graph.V 5 >>> cobalt_complex = io.experimental.from_smiles("Br[Co@OH12](Cl)(I)(F)(S)C") >>> cobalt_index = cobalt_complex.graph.atoms_of_element(utils.ElementType.Co)[0] >>> permutator = cobalt_complex.stereopermutators.option(cobalt_index) >>> permutator is not None True >>> permutator.assigned is not None True
-
scine_molassembler.io.experimental.
from_smiles_multiple
(smiles_str: str) → List[scine_molassembler.Molecule]¶ Parse a smiles string containing possibly multiple molecules
The smiles parser is implemented according to the OpenSMILES spec. It supports the following features:
Arbitrarily many molecules in a string
Isotope markers
Valence filling of the organic subset
Set shapes from VSEPR
Ring closures
Stereo markers - Double bond - Tetrahedral - Square planar - Trigonal bipyramidal - Octahedral
- Parameters
smiles_str – A smiles string containing possibly multiple molecules
>>> methane_and_ammonia = io.experimental.from_smiles_multiple("C.[NH4+]") >>> len(methane_and_ammonia) == 2 True