openforcefield.topology.
Molecule
(*args, **kwargs)[source]¶Mutable chemical representation of a molecule, such as a small molecule or biopolymer.
Examples
Create a molecule from an sdf file
>>> from openforcefield.utils import get_data_filename
>>> sdf_filepath = get_data_filename('molecules/ethanol.sdf')
>>> molecule = Molecule(sdf_filepath)
Convert to OpenEye OEMol object
>>> oemol = molecule.to_openeye()
Create a molecule from an OpenEye molecule
>>> molecule = Molecule.from_openeye(oemol)
Convert to RDKit Mol object
>>> rdmol = molecule.to_rdkit()
Create a molecule from an RDKit molecule
>>> molecule = Molecule.from_rdkit(rdmol)
Create a molecule from IUPAC name (requires the OpenEye toolkit)
>>> molecule = Molecule.from_iupac('imatinib')
Create a molecule from SMILES
>>> molecule = Molecule.from_smiles('Cc1ccccc1')
Warning
This API is experimental and subject to change.
Attributes: |
|
---|
Methods
add_atom (atomic_number, formal_charge, …) |
Add an atom |
add_bond (atom1, atom2, bond_order, is_aromatic) |
Add a bond between two specified atom indices |
add_bond_charge_virtual_site (atoms, distance) |
Create a bond charge-type virtual site, in which the location of the charge is specified by the position of two atoms. |
add_conformer (coordinates) |
# TODO: Should this not be public? Adds a conformer of the molecule |
add_divalent_lone_pair_virtual_site (atoms, …) |
Create a divalent lone pair-type virtual site, in which the location of the charge is specified by the position of three atoms. |
add_monovalent_lone_pair_virtual_site (atoms, …) |
Create a bond charge-type virtual site, in which the location of the charge is specified by the position of three atoms. |
add_trivalent_lone_pair_virtual_site (atoms, …) |
Create a trivalent lone pair-type virtual site, in which the location of the charge is specified by the position of four atoms. |
chemical_environment_matches (query[, …]) |
Retrieve all matches for a given chemical environment query. |
compute_partial_charges ([toolkit_registry]) |
Warning! Not Implemented! Calculate partial atomic charges for this molecule using an underlying toolkit |
compute_partial_charges_am1bcc ([…]) |
Calculate partial atomic charges for this molecule using AM1-BCC run by an underlying toolkit |
compute_wiberg_bond_orders ([charge_model, …]) |
Calculate wiberg bond orders for this molecule using an underlying toolkit |
from_bson (serialized) |
Instantiate an object from a BSON serialized representation. |
from_dict (molecule_dict) |
Create a new Molecule from a dictionary representation |
from_file (filename[, file_format, …]) |
Create one or more molecules from a file |
from_iupac (iupac_name, **kwargs) |
Generate a molecule from IUPAC or common name |
from_json (serialized) |
Instantiate an object from a JSON serialized representation. |
from_messagepack (serialized) |
Instantiate an object from a MessagePack serialized representation. |
from_openeye (oemol[, allow_undefined_stereo]) |
Create a Molecule from an OpenEye molecule. |
from_pickle (serialized) |
Instantiate an object from a pickle serialized representation. |
from_rdkit (rdmol[, allow_undefined_stereo]) |
Create a Molecule from an RDKit molecule. |
from_smiles (smiles[, …]) |
Construct a Molecule from a SMILES representation |
from_toml (serialized) |
Instantiate an object from a TOML serialized representation. |
from_topology (topology) |
Return a Molecule representation of an openforcefield Topology containing a single Molecule object. |
from_xml (serialized) |
Instantiate an object from an XML serialized representation. |
from_yaml (serialized) |
Instantiate from a YAML serialized representation. |
generate_conformers ([toolkit_registry, …]) |
Generate conformers for this molecule using an underlying toolkit |
get_fractional_bond_orders ([method, …]) |
Get fractional bond orders. |
is_isomorphic (other[, …]) |
Determines whether the molecules are isomorphic by comparing their graphs. |
to_bson () |
Return a BSON serialized representation. |
to_dict () |
Return a dictionary representation of the molecule. |
to_file (outfile, outfile_format[, …]) |
Write the current molecule to a file or file-like object |
to_iupac () |
Generate IUPAC name from Molecule |
to_json ([indent]) |
Return a JSON serialized representation. |
to_messagepack () |
Return a MessagePack representation. |
to_networkx () |
Generate a NetworkX undirected graph from the Molecule. |
to_openeye ([aromaticity_model]) |
Create an OpenEye molecule |
to_pickle () |
Return a pickle serialized representation. |
to_rdkit ([aromaticity_model]) |
Create an RDKit molecule |
to_smiles ([toolkit_registry]) |
Return a canonical isomeric SMILES representation of the current molecule |
to_toml () |
Return a TOML serialized representation. |
to_topology () |
Return an openforcefield Topology representation containing one copy of this molecule |
to_xml ([indent]) |
Return an XML representation. |
to_yaml () |
Return a YAML serialized representation. |
__init__
(*args, **kwargs)[source]¶Create a new Molecule object
Parameters: |
|
---|
Examples
Create an empty molecule:
>>> empty_molecule = Molecule()
Create a molecule from a file that can be used to construct a molecule, using either a filename or file-like object:
>>> from openforcefield.utils import get_data_filename
>>> sdf_filepath = get_data_filename('molecules/ethanol.sdf')
>>> molecule = Molecule(sdf_filepath)
>>> molecule = Molecule(open(sdf_filepath, 'r'), file_format='sdf')
>>> import gzip
>>> mol2_gz_filepath = get_data_filename('molecules/toluene.mol2.gz')
>>> molecule = Molecule(gzip.GzipFile(mol2_gz_filepath, 'r'), file_format='mol2')
Create a molecule from another molecule:
>>> molecule_copy = Molecule(molecule)
Convert to OpenEye OEMol object
>>> oemol = molecule.to_openeye()
Create a molecule from an OpenEye molecule:
>>> molecule = Molecule(oemol)
Convert to RDKit Mol object
>>> rdmol = molecule.to_rdkit()
Create a molecule from an RDKit molecule:
>>> molecule = Molecule(rdmol)
Create a molecule from a serialized molecule object:
>>> serialized_molecule = molecule.__getstate__()
>>> molecule_copy = Molecule(serialized_molecule)
Methods
__init__ (*args, **kwargs) |
Create a new Molecule object |
add_atom (atomic_number, formal_charge, …) |
Add an atom |
add_bond (atom1, atom2, bond_order, is_aromatic) |
Add a bond between two specified atom indices |
add_bond_charge_virtual_site (atoms, distance) |
Create a bond charge-type virtual site, in which the location of the charge is specified by the position of two atoms. |
add_conformer (coordinates) |
# TODO: Should this not be public? Adds a conformer of the molecule |
add_divalent_lone_pair_virtual_site (atoms, …) |
Create a divalent lone pair-type virtual site, in which the location of the charge is specified by the position of three atoms. |
add_monovalent_lone_pair_virtual_site (atoms, …) |
Create a bond charge-type virtual site, in which the location of the charge is specified by the position of three atoms. |
add_trivalent_lone_pair_virtual_site (atoms, …) |
Create a trivalent lone pair-type virtual site, in which the location of the charge is specified by the position of four atoms. |
chemical_environment_matches (query[, …]) |
Retrieve all matches for a given chemical environment query. |
compute_partial_charges ([toolkit_registry]) |
Warning! Not Implemented! Calculate partial atomic charges for this molecule using an underlying toolkit |
compute_partial_charges_am1bcc ([…]) |
Calculate partial atomic charges for this molecule using AM1-BCC run by an underlying toolkit |
compute_wiberg_bond_orders ([charge_model, …]) |
Calculate wiberg bond orders for this molecule using an underlying toolkit |
from_bson (serialized) |
Instantiate an object from a BSON serialized representation. |
from_dict (molecule_dict) |
Create a new Molecule from a dictionary representation |
from_file (filename[, file_format, …]) |
Create one or more molecules from a file |
from_iupac (iupac_name, **kwargs) |
Generate a molecule from IUPAC or common name |
from_json (serialized) |
Instantiate an object from a JSON serialized representation. |
from_messagepack (serialized) |
Instantiate an object from a MessagePack serialized representation. |
from_openeye (oemol[, allow_undefined_stereo]) |
Create a Molecule from an OpenEye molecule. |
from_pickle (serialized) |
Instantiate an object from a pickle serialized representation. |
from_rdkit (rdmol[, allow_undefined_stereo]) |
Create a Molecule from an RDKit molecule. |
from_smiles (smiles[, …]) |
Construct a Molecule from a SMILES representation |
from_toml (serialized) |
Instantiate an object from a TOML serialized representation. |
from_topology (topology) |
Return a Molecule representation of an openforcefield Topology containing a single Molecule object. |
from_xml (serialized) |
Instantiate an object from an XML serialized representation. |
from_yaml (serialized) |
Instantiate from a YAML serialized representation. |
generate_conformers ([toolkit_registry, …]) |
Generate conformers for this molecule using an underlying toolkit |
get_fractional_bond_orders ([method, …]) |
Get fractional bond orders. |
is_isomorphic (other[, …]) |
Determines whether the molecules are isomorphic by comparing their graphs. |
to_bson () |
Return a BSON serialized representation. |
to_dict () |
Return a dictionary representation of the molecule. |
to_file (outfile, outfile_format[, …]) |
Write the current molecule to a file or file-like object |
to_iupac () |
Generate IUPAC name from Molecule |
to_json ([indent]) |
Return a JSON serialized representation. |
to_messagepack () |
Return a MessagePack representation. |
to_networkx () |
Generate a NetworkX undirected graph from the Molecule. |
to_openeye ([aromaticity_model]) |
Create an OpenEye molecule |
to_pickle () |
Return a pickle serialized representation. |
to_rdkit ([aromaticity_model]) |
Create an RDKit molecule |
to_smiles ([toolkit_registry]) |
Return a canonical isomeric SMILES representation of the current molecule |
to_toml () |
Return a TOML serialized representation. |
to_topology () |
Return an openforcefield Topology representation containing one copy of this molecule |
to_xml ([indent]) |
Return an XML representation. |
to_yaml () |
Return a YAML serialized representation. |
Attributes
angles |
Get an iterator over all i-j-k angles. |
atoms |
Iterate over all Atom objects. |
bonds |
Iterate over all Bond objects. |
conformers |
Iterate over all conformers in this molecule. |
impropers |
Iterate over all proper torsions in the molecule |
n_angles |
int: number of angles in the Molecule. |
n_atoms |
The number of Atom objects. |
n_bonds |
The number of Bond objects. |
n_conformers |
Iterate over all Atom objects. |
n_impropers |
int: number of improper torsions in the Molecule. |
n_particles |
The number of Particle objects, which corresponds to how many positions must be used. |
n_propers |
int: number of proper torsions in the Molecule. |
n_virtual_sites |
The number of VirtualSite objects. |
name |
The name (or title) of the molecule |
partial_charges |
Returns the partial charges (if present) on the molecule |
particles |
Iterate over all Particle objects. |
propers |
Iterate over all proper torsions in the molecule |
properties |
The properties dictionary of the molecule |
torsions |
Get an iterator over all i-j-k-l torsions. |
total_charge |
Return the total charge on the molecule |
virtual_sites |
Iterate over all VirtualSite objects. |
angles
¶Get an iterator over all i-j-k angles.
atoms
¶Iterate over all Atom objects.
bonds
¶Iterate over all Bond objects.
chemical_environment_matches
(query, toolkit_registry=<openforcefield.utils.toolkits.ToolkitRegistry object>)¶Retrieve all matches for a given chemical environment query.
Parameters: |
|
---|---|
Returns: |
|
Examples
Retrieve all the carbon-carbon bond matches in a molecule
>>> molecule = Molecule.from_iupac('imatinib')
>>> matches = molecule.chemical_environment_matches('[#6X3:1]~[#6X3:2]')
compute_partial_charges
(toolkit_registry=<openforcefield.utils.toolkits.ToolkitRegistry object>)¶Warning! Not Implemented! Calculate partial atomic charges for this molecule using an underlying toolkit
Parameters: |
|
---|---|
Raises: |
|
Examples
>>> molecule = Molecule.from_smiles('CCCCCC')
>>> molecule.generate_conformers()
>>> molecule.compute_partial_charges()
compute_partial_charges_am1bcc
(toolkit_registry=<openforcefield.utils.toolkits.ToolkitRegistry object>)¶Calculate partial atomic charges for this molecule using AM1-BCC run by an underlying toolkit
Parameters: |
|
---|---|
Raises: |
|
Examples
>>> molecule = Molecule.from_smiles('CCCCCC')
>>> molecule.generate_conformers()
>>> molecule.compute_partial_charges_am1bcc()
compute_wiberg_bond_orders
(charge_model=None, toolkit_registry=<openforcefield.utils.toolkits.ToolkitRegistry object>)¶Calculate wiberg bond orders for this molecule using an underlying toolkit
Parameters: |
|
---|---|
Raises: |
|
Examples
>>> molecule = Molecule.from_smiles('CCCCCC')
>>> molecule.generate_conformers()
>>> molecule.compute_wiberg_bond_orders()
conformers
¶Iterate over all conformers in this molecule.
from_bson
(serialized)¶Instantiate an object from a BSON serialized representation.
Specification: http://bsonspec.org/
Parameters: |
|
---|---|
Returns: |
|
from_dict
(molecule_dict)¶Create a new Molecule from a dictionary representation
Parameters: |
|
---|---|
Returns: |
|
from_file
(filename, file_format=None, toolkit_registry=<openforcefield.utils.toolkits.ToolkitRegistry object>, allow_undefined_stereo=False)¶Create one or more molecules from a file
Parameters: |
|
---|---|
Returns: |
|
Examples
>>> from openforcefield.tests.utils import get_monomer_mol2file
>>> mol2_filename = get_monomer_mol2file('cyclohexane')
>>> molecule = Molecule.from_file(mol2_filename)
from_iupac
(iupac_name, **kwargs)¶Generate a molecule from IUPAC or common name
Parameters: |
|
---|---|
Returns: |
|
Examples
Create a molecule from a common name
>>> molecule = Molecule.from_iupac('4-[(4-methylpiperazin-1-yl)methyl]-N-(4-methyl-3-{[4-(pyridin-3-yl)pyrimidin-2-yl]amino}phenyl)benzamide')
Create a molecule from a common name
>>> molecule = Molecule.from_iupac('imatinib')
from_json
(serialized)¶Instantiate an object from a JSON serialized representation.
Specification: https://www.json.org/
Parameters: |
|
---|---|
Returns: |
|
from_messagepack
(serialized)¶Instantiate an object from a MessagePack serialized representation.
Specification: https://msgpack.org/index.html
Parameters: |
|
---|---|
Returns: |
|
from_openeye
(oemol, allow_undefined_stereo=False)¶Create a Molecule from an OpenEye molecule.
Requires the OpenEye toolkit to be installed.
Parameters: |
|
---|---|
Returns: |
|
Examples
Create a Molecule from an OpenEye OEMol
>>> from openeye import oechem
>>> from openforcefield.tests.utils import get_data_filename
>>> ifs = oechem.oemolistream(get_data_filename('systems/monomers/ethanol.mol2'))
>>> oemols = list(ifs.GetOEGraphMols())
>>> molecule = Molecule.from_openeye(oemols[0])
from_pickle
(serialized)¶Instantiate an object from a pickle serialized representation.
Warning
This is not recommended for safe, stable storage since the pickle specification may change between Python versions.
Parameters: |
|
---|---|
Returns: |
|
from_rdkit
(rdmol, allow_undefined_stereo=False)¶Create a Molecule from an RDKit molecule.
Requires the RDKit to be installed.
Parameters: |
|
---|---|
Returns: |
|
Examples
Create a molecule from an RDKit molecule
>>> from rdkit import Chem
>>> from openforcefield.tests.utils import get_data_filename
>>> rdmol = Chem.MolFromMolFile(get_data_filename('systems/monomers/ethanol.sdf'))
>>> molecule = Molecule.from_rdkit(rdmol)
from_smiles
(smiles, hydrogens_are_explicit=False, toolkit_registry=<openforcefield.utils.toolkits.ToolkitRegistry object>)¶Construct a Molecule from a SMILES representation
Parameters: |
|
---|---|
Returns: |
|
Examples
>>> molecule = Molecule.from_smiles('Cc1ccccc1')
from_toml
(serialized)¶Instantiate an object from a TOML serialized representation.
Specification: https://github.com/toml-lang/toml
Parameters: |
|
---|---|
Returns: |
|
from_topology
(topology)¶Return a Molecule representation of an openforcefield Topology containing a single Molecule object.
Parameters: | |
---|---|
Returns: |
|
Raises: |
|
Examples
Create a molecule from a Topology object that contains exactly one molecule
>>> molecule = Molecule.from_topology(topology) # doctest: +SKIP
from_xml
(serialized)¶Instantiate an object from an XML serialized representation.
Specification: https://www.w3.org/XML/
Parameters: |
|
---|---|
Returns: |
|
from_yaml
(serialized)¶Instantiate from a YAML serialized representation.
Specification: http://yaml.org/
Parameters: |
|
---|---|
Returns: |
|
generate_conformers
(toolkit_registry=<openforcefield.utils.toolkits.ToolkitRegistry object>, num_conformers=10, clear_existing=True)¶Generate conformers for this molecule using an underlying toolkit
Parameters: |
|
---|---|
Raises: |
|
Examples
>>> molecule = Molecule.from_smiles('CCCCCC')
>>> molecule.generate_conformers()
get_fractional_bond_orders
(method='Wiberg', toolkit_registry=<openforcefield.utils.toolkits.ToolkitRegistry object>)¶Get fractional bond orders.
Examples
Get fractional Wiberg bond orders
>>> molecule = Molecule.from_iupac('imatinib')
>>> molecule.generate_conformers()
>>> fractional_bond_orders = molecule.get_fractional_bond_orders(method='Wiberg')
impropers
¶Iterate over all proper torsions in the molecule
is_isomorphic
(other, compare_atom_stereochemistry=True, compare_bond_stereochemistry=True)¶Determines whether the molecules are isomorphic by comparing their graphs.
Parameters: |
|
---|---|
Returns: |
|
n_angles
¶int: number of angles in the Molecule.
n_atoms
¶The number of Atom objects.
n_bonds
¶The number of Bond objects.
n_conformers
¶Iterate over all Atom objects.
n_impropers
¶int: number of improper torsions in the Molecule.
n_particles
¶The number of Particle objects, which corresponds to how many positions must be used.
n_propers
¶int: number of proper torsions in the Molecule.
n_virtual_sites
¶The number of VirtualSite objects.
name
¶The name (or title) of the molecule
partial_charges
¶Returns the partial charges (if present) on the molecule
Returns: |
|
---|
particles
¶Iterate over all Particle objects.
propers
¶Iterate over all proper torsions in the molecule
properties
¶The properties dictionary of the molecule
to_bson
()¶Return a BSON serialized representation.
Specification: http://bsonspec.org/
Returns: |
|
---|
to_dict
()¶Return a dictionary representation of the molecule.
Returns: |
|
---|
to_file
(outfile, outfile_format, toolkit_registry=<openforcefield.utils.toolkits.ToolkitRegistry object>)¶Write the current molecule to a file or file-like object
Parameters: |
|
---|---|
Raises: |
|
Examples
>>> molecule = Molecule.from_iupac('imatinib')
>>> molecule.to_file('imatinib.mol2', outfile_format='mol2') # doctest: +SKIP
>>> molecule.to_file('imatinib.sdf', outfile_format='sdf') # doctest: +SKIP
>>> molecule.to_file('imatinib.pdb', outfile_format='pdb') # doctest: +SKIP
to_iupac
()¶Generate IUPAC name from Molecule
Returns: |
|
---|
Examples
>>> from openforcefield.utils import get_data_filename
>>> sdf_filepath = get_data_filename('molecules/ethanol.sdf')
>>> molecule = Molecule(sdf_filepath)
>>> iupac_name = molecule.to_iupac()
to_json
(indent=None)¶Return a JSON serialized representation.
Specification: https://www.json.org/
Parameters: |
|
---|---|
Returns: |
|
to_messagepack
()¶Return a MessagePack representation.
Specification: https://msgpack.org/index.html
Returns: |
|
---|
to_networkx
()¶Generate a NetworkX undirected graph from the Molecule.
Nodes are Atoms labeled with particle indices and atomic elements (via the element
node atrribute).
Edges denote chemical bonds between Atoms.
Virtual sites are not included, since they lack a concept of chemical connectivity.
Returns: |
|
---|
Examples
Retrieve the bond graph for imatinib (OpenEye toolkit required)
>>> molecule = Molecule.from_iupac('imatinib')
>>> nxgraph = molecule.to_networkx()
to_openeye
(aromaticity_model='OEAroModel_MDL')¶Create an OpenEye molecule
Requires the OpenEye toolkit to be installed.
Parameters: |
|
---|---|
Returns: |
|
Examples
Create an OpenEye molecule from a Molecule
>>> molecule = Molecule.from_smiles('CC')
>>> oemol = molecule.to_openeye()
to_pickle
()¶Return a pickle serialized representation.
Warning
This is not recommended for safe, stable storage since the pickle specification may change between Python versions.
Returns: |
|
---|
to_rdkit
(aromaticity_model='OEAroModel_MDL')¶Create an RDKit molecule
Requires the RDKit to be installed.
Parameters: |
|
---|---|
Returns: |
|
Examples
Convert a molecule to RDKit
>>> from openforcefield.utils import get_data_filename
>>> sdf_filepath = get_data_filename('molecules/ethanol.sdf')
>>> molecule = Molecule(sdf_filepath)
>>> rdmol = molecule.to_rdkit()
to_smiles
(toolkit_registry=<openforcefield.utils.toolkits.ToolkitRegistry object>)¶Return a canonical isomeric SMILES representation of the current molecule
Note
RDKit and OpenEye versions will not necessarily return the same representation.
Parameters: |
|
---|---|
Returns: |
|
Examples
>>> from openforcefield.utils import get_data_filename
>>> sdf_filepath = get_data_filename('molecules/ethanol.sdf')
>>> molecule = Molecule(sdf_filepath)
>>> smiles = molecule.to_smiles()
to_toml
()¶Return a TOML serialized representation.
Specification: https://github.com/toml-lang/toml
Returns: |
|
---|
to_topology
()¶Return an openforcefield Topology representation containing one copy of this molecule
Returns: |
|
---|
Examples
>>> molecule = Molecule.from_iupac('imatinib')
>>> topology = molecule.to_topology()
to_xml
(indent=2)¶Return an XML representation.
Specification: https://www.w3.org/XML/
Parameters: |
|
---|---|
Returns: |
|
to_yaml
()¶Return a YAML serialized representation.
Specification: http://yaml.org/
Returns: |
|
---|
torsions
¶Get an iterator over all i-j-k-l torsions. Note that i-j-k-i torsions (cycles) are excluded.
Returns: |
|
---|
total_charge
¶Return the total charge on the molecule
virtual_sites
¶Iterate over all VirtualSite objects.
add_atom
(atomic_number, formal_charge, is_aromatic, stereochemistry=None, name=None)[source]¶Add an atom
Parameters: |
|
---|---|
Returns: |
|
Examples
Define a methane molecule
>>> molecule = Molecule()
>>> molecule.name = 'methane'
>>> C = molecule.add_atom(6, 0, False)
>>> H1 = molecule.add_atom(1, 0, False)
>>> H2 = molecule.add_atom(1, 0, False)
>>> H3 = molecule.add_atom(1, 0, False)
>>> H4 = molecule.add_atom(1, 0, False)
>>> bond_idx = molecule.add_bond(C, H1, False, 1)
>>> bond_idx = molecule.add_bond(C, H2, False, 1)
>>> bond_idx = molecule.add_bond(C, H3, False, 1)
>>> bond_idx = molecule.add_bond(C, H4, False, 1)
add_bond_charge_virtual_site
(atoms, distance, charge_increments=None, weights=None, epsilon=None, sigma=None, rmin_half=None, name='')[source]¶Create a bond charge-type virtual site, in which the location of the charge is specified by the position of two atoms. This supports placement of a virtual site S along a vector between two specified atoms, e.g. to allow for a sigma hole for halogens or similar contexts. With positive values of the distance, the virtual site lies outside the first indexed atom. Parameters ———- atoms : list of openforcefield.topology.molecule.Atom objects or ints of shape [N
The atoms defining the virtual site’s position or their indices
distance : float
Returns: |
|
---|
add_monovalent_lone_pair_virtual_site
(atoms, distance, out_of_plane_angle, in_plane_angle, **kwargs)[source]¶Create a bond charge-type virtual site, in which the location of the charge is specified by the position of three atoms.
TODO : Do “weights” have any meaning here?
Parameters: |
|
---|---|
Returns: |
|
add_divalent_lone_pair_virtual_site
(atoms, distance, out_of_plane_angle, in_plane_angle, **kwargs)[source]¶Create a divalent lone pair-type virtual site, in which the location of the charge is specified by the position of three atoms.
TODO : Do “weights” have any meaning here?
Parameters: |
|
---|---|
Returns: |
|
add_trivalent_lone_pair_virtual_site
(atoms, distance, out_of_plane_angle, in_plane_angle, **kwargs)[source]¶Create a trivalent lone pair-type virtual site, in which the location of the charge is specified by the position of four atoms.
TODO : Do “weights” have any meaning here?
Parameters: |
|
---|---|
Returns: |
|
add_bond
(atom1, atom2, bond_order, is_aromatic, stereochemistry=None, fractional_bond_order=None)[source]¶Add a bond between two specified atom indices
Parameters: |
|
---|---|
Returns: |
|
add_conformer
(coordinates)[source]¶# TODO: Should this not be public? Adds a conformer of the molecule
Parameters: |
|
---|