morfeus.conformer module#
Conformer tools.
- class morfeus.conformer.Conformer(elements, coordinates, energy=None, degeneracy=1, properties=None, cip_label=None)[source]#
Bases:
object
Conformer with coordinates, energy and properties.
- Parameters:
elements (Iterable[int] | Iterable[str]) – Elements as atomic symbols or numbers
coordinates (ArrayLike2D) – Coordinates (Å)
energy (float | None) – Energy (a.u.)
degeneracy (int) – Degeneracy
properties (dict[str, float] | None) – Conformers properties.
cip_label (tuple[str] | None) – Tuple of CIP labels for all atoms: ‘R’, ‘S’ or ‘’
- cip_label#
Tuple of CIP labels for all atoms: ‘R’, ‘S’ or ‘’
- Type:
tuple[str, …] | None
- coordinates#
Coordinates (Å)
- Type:
Array2DFloat
- degeneracy#
Degeneracy
- Type:
int
- elements#
Elements as atomic symbols or numbers
- Type:
Array1DFloat
- energy#
Energy (a.u.)
- Type:
float | None
- properties#
Conformers properties.
- Type:
dict[str, float]
- class morfeus.conformer.ConformerEnsemble(elements, conformer_coordinates=None, energies=None, connectivity_matrix=None, degeneracies=None, properties=None, charge=None, multiplicity=None, formal_charges=None, ref_cip_label=None)[source]#
Bases:
object
Conformer ensemble object.
Supports sorting, pruning, optimization and single-point calculations.
- Parameters:
elements (Iterable[int] | Iterable[str]) – Elements as atomic symbols or numbers
conformer_coordinates (ArrayLike3D | None) – Conformer coordinates (Å)
energies (ArrayLike1D | None) – Energy (a.u.)
connectivity_matrix (ArrayLike2D | None) – Connectivity matrix
degeneracies (ArrayLike1D | None) – Degeneracies
properties (Mapping[str, ArrayLike1D] | None) – Conformers properties.
charge (int | None) – Molecular charge
multiplicity (int | None) – Molecular multiplicity.
formal_charges (ArrayLike1D | None) – Atomic formal charges
ref_cip_label (tuple[str, ...] | None) – Tuple of CIP labels for all atoms: ‘R’, ‘S’ or ‘’
- charge#
Molecular charge
- Type:
int
- conformers#
Conformers
- connectivity_matrix#
Connectivity matrix
- Type:
Array2DInt
- elements#
Elements as atomic symbols or numbers
- Type:
Array1DInt
- formal_charges#
Atomic formal charges.
- Type:
Array1DInt
- mol#
RDKit mol object.
- Type:
Chem.Mol
- multiplicity#
Molecular multiplicity.
- Type:
int
- ref_cip_label#
Tuple of CIP labels for all atoms: ‘R’, ‘S’ or ‘’
- Type:
tuple[str, …] | None
- add_conformers(conformer_coordinates, energies=None, degeneracies=None, properties=None)[source]#
Add conformer to ensemble.
- Parameters:
conformer_coordinates (ArrayLike3D) – Conformer coordinates (Å)
energies (ArrayLike1D | None) – Energies (a.u.)
degeneracies (ArrayLike1D | None) – Degeneracies
properties (Mapping[str, ArrayLike1D] | None) – Conformer properties.
- Returns:
Self
- Return type:
self
- add_inverted()[source]#
Add inverted images of all conformers.
Scrambles stereochemistry and leads to redundant conformers so use with care and prune based on RMSD as postprocessing.
- Returns:
Self
- Return type:
self
- boltzmann_average_dT(property_name, temperature=298.15)[source]#
Calculate temperature derivative of Boltzmann average of property.
- Parameters:
property_name (str) – Name of property
temperature (float) – Temperature (K)
- Returns:
Derivative of Boltzmann average
- Return type:
derivative
- boltzmann_statistic(property_name, temperature=298.15, statistic='avg')[source]#
Calculate Boltzmann staistic of property over ensemble.
- Parameters:
property_name (str) – Name of property
statistic (str) – Boltzmann statistic: ‘avg’, ‘var’ or ‘std’
temperature (float) – Temperature (K)
- Returns:
Boltzmann statistic
- Return type:
statistic
- boltzmann_weights(temperature=298.15)[source]#
Calculate Boltzmann weights for ensemble.
- Parameters:
temperature (float) – Temperature (K)
- Returns:
Conformer weights (normalized to unity)
- Return type:
weights
- condense_enantiomeric(thres=None)[source]#
Condense enantiomers into single enantiomer per pair.
- Parameters:
thres (float | None) – RMSD threshold for assessing enantiomers.
- Returns:
Self
- Return type:
self
- Raises:
Exception – If molecule is chiral.
- detect_enantiomers(thres=0.01, method='rdkit', include_hs=False)[source]#
Detect enantiomers in ensemble.
- Parameters:
thres (float) – RMSD threshold for detecting enantiomers in terms of coordinates
method (str) – RMSD calculation method: ‘obrms-batch’, ‘obrms-iter’, openbabel’ or ‘spyrmsd’
include_hs (bool) – Whether to include H atoms when determining enantiomers
- Returns:
- Mapping of enantiomer with conformer id as keys and enantiomer
ids as values.
- Return type:
enantiomers
- classmethod from_crest(path)[source]#
Generate conformer ensemble from CREST output.
- Parameters:
path (str | PathLike) – Path to CREST folder
- Returns:
Conformer ensemble object.
- Return type:
ce
- classmethod from_ob_ff(*args, generate_rdkit_mol=False, update_charges=True, update_connectivity=True, update_multiplicity=True, **kwargs)[source]#
Generate conformer ensemble from OpenBabel with FF method.
See the documentation for the function conformers_from_ob_ff for more information.
- Parameters:
*args (Any) – Positional arguments for conformers_from_ob_ga
generate_rdkit_mol (bool) – Generate RDKit mol object for ensemble
update_charges (bool) – Update formal charges from generated RDKit Mol object. Only used if generate_rdkit_mol is True.
update_connectivity (bool) – Update connectivity from generated RDKit Mol object. Only used if generate_rdkit_mol is True.
update_multiplicity (bool) – Update multiplicity from generated RDKit Mol object. Only used if generate_rdkit_mol is True.
**kwargs (Any) – Keyword arguments for conformers_from_ob_ga
- Returns:
Conformer ensemble.
- Return type:
ce
- classmethod from_ob_ga(*args, generate_rdkit_mol=False, update_charges=True, update_connectivity=True, update_multiplicity=True, **kwargs)[source]#
Generate conformer ensemble from OpenBabel with GA method.
See the documentation for the function conformers_from_ob_ga for more information.
- Parameters:
*args (Any) – Positional arguments for conformers_from_ob_ga
generate_rdkit_mol (bool) – Generate RDKit mol object for ensemble
update_charges (bool) – Update formal charges from generated RDKit Mol object. Only used if generate_rdkit_mol is True.
update_connectivity (bool) – Update connectivity from generated RDKit Mol object. Only used if generate_rdkit_mol is True.
update_multiplicity (bool) – Update multiplicity from generated RDKit Mol object. Only used if generate_rdkit_mol is True.
**kwargs (Any) – Keyword arguments for conformers_from_ob_ga
- Returns:
Conformer ensemble.
- Return type:
ce
- classmethod from_rdkit(*args, **kwargs)[source]#
Generate conformer ensemble from RDKit.
See the documentation for the function conformers_from_rdkit for more information.
- Parameters:
*args (Any) – Positional arguments
**kwargs (Any) – Keyword arguments
- Returns:
Conformer ensemble object.
- Return type:
ce
- generate_mol(update_charges=True, update_connectivity=True, update_multiplicity=True)[source]#
Generate RDKit Mol object.
- Parameters:
update_charges (bool) – Update formal charges from generated RDKit Mol object
update_connectivity (bool) – Update connectivity from generated RDKit Mol object
update_multiplicity (bool) – Update multiplicity from generated RDKit Mol object
- Return type:
None
- get_cip_labels()[source]#
Generate tuples of CIP labels for conformer.
- Returns:
Tuples of CIP labels for each conformer.
- Return type:
cip_labels
- get_coordinates()[source]#
Get conformer coordinates.
- Returns:
Conformer coordinates (Å)
- Return type:
conformer_coordinates
- get_degeneracies()[source]#
Get conformer degeneracies.
- Returns:
Degeneracies
- Return type:
degeneriacies
- get_properties()[source]#
Get conformer properties.
- Returns:
Conformer properties
- Return type:
properties
- get_relative_energies(unit='kcal/mol', relative=True)[source]#
Get conformer energies with choice of units and reference value.
- Parameters:
unit (str) – Unit of returned energies: ‘hartree’, ‘kcal/mol’ or ‘kJ/mol’.
relative (bool) – Return energies relative to lowest-energy conformer.
- Returns:
Conformer energies.
- Return type:
energies
- get_rmsd(i_s=None, j_s=None, include_hs=False, symmetry=False, method='rdkit')[source]#
Get RSMD between conformers.
For very small systems ‘openbabel’ or ‘spyrmsd’ work well. For larger systems a significant speed-up is attained with ‘rdkit’, ‘obrms-batch’ or ‘obrms-iter’.
- Parameters:
i_s (ArrayLike1D | None) – Indices of conformers
j_s (ArrayLike1D | None) – Indices of conformers
include_hs (bool) – Whether to include H atoms in RMSD calculation. Ignored for ‘obrms-iter’ and ‘obrms-batch’ that only use heavy atoms.
symmetry (bool) – Consider symmetry (requires connectivity matrix). Ignored for ‘obrms-iter’ and ‘obrms-batch’ that always use symmetry.
method (str) – RMSD calculation method: ‘obrms-batch’, ‘obrms-iter’, ‘openbabel’, ‘rdkit’ or ‘spyrmsd’.
- Returns:
RSMDs (Å)
- Return type:
rmsds
- Raises:
ValueError – If method not supported.
- property n_conformers: int#
Number of conformers.
- optimize_qc_engine(ids=None, program=None, model=None, keywords=None, local_options=None, procedure='berny')[source]#
Optimize conformers with QCEngine interface.
- Parameters:
ids (Sequence[int] | None) – Conformer indices (1-indexed). If None, all are optimized
program (str | None) – QCEngine program
model (dict[str, Any] | None) – QCEngine model
keywords (dict[str, Any] | None) – QCEngine keywords
local_options (dict[str, Any] | None) – QCEngine local options
procedure (str) – QCEngine procedure
- Returns:
Self
- Return type:
self
- Raises:
Exception – When RDKit requested without formal charges.
- prune_enantiomers(keep='original', ref_label=None)[source]#
Prune conformers so that only one enantiomer is present in the ensemble.
- Parameters:
keep (str) – Which enantiomer to keep: ‘original’, ‘most common’ or ‘specified’. Choice of ‘original’ requires that the ref_cip_label attribute is set. Choice of ‘specified’ requires ref_label to be given.
ref_label (tuple[str, ...] | None) – Reference CIP labels for all atoms
- Returns:
Self
- Return type:
self
- prune_energy(threshold=3.0, unit='kcal/mol')[source]#
Prune conformers based on energy compared to minimum energy conformer.
- Parameters:
threshold (float) – Energy threshold for pruning
unit (str) – Unit for energy threshold ‘hartree’, ‘kcal/mol’ or ‘kJ/mol’
- Returns:
Self
- Return type:
self
- prune_rmsd(thres=0.35, include_hs=False, symmetry=False, method='rdkit')[source]#
Prune conformers based on RMSD.
- Parameters:
thres (float) – Threshold for RSMD pruning (Å)
include_hs (bool) – Whether to include H atoms in RMSD calculation. Ignored for ‘obrms-iter’ and ‘obrms-batch’ which only use heavy atoms.
symmetry (bool) – Consider symmetry (requires connectivity matrix). Ignored for ‘obrms-iter’ and ‘obrms-batch’ which always use symmetry.
method (str) – RMSD calculation method: ‘obrms-batch’, ‘obrms-iter’, ‘openbabel’, ‘rdkit’ or ‘spyrmsd’
- Returns:
Self
- Return type:
self
- set_coordinates(conformer_coordinates)[source]#
Set conformer coordinates.
- Parameters:
conformer_coordinates (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – Conformer coordinates (Å)
- Returns:
Self
- Return type:
self
- Raises:
ValueError – When Number of conformer coordinates is different from number of conformers
- set_degeneracies(degeneracies)[source]#
Set degeneracies.
- Parameters:
degeneracies (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – Degeneracies
- Returns:
Self
- Return type:
self
- Raises:
ValueError – When number of degeneracies is different from number of conformers
- set_energies(energies)[source]#
Set energies.
- Parameters:
energies (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – Energy (a.u.)
- Returns:
Self
- Return type:
self
- Raises:
ValueError – When number of energies is different from number of conformers
- set_multiplicity_from_mol()[source]#
Sets multiplicity based on unpaired electrons in Mol object.
- Return type:
- set_properties(key, values)[source]#
Set conformer properties.
- Parameters:
key (str) – Name of property
values (Iterable[float]) – Property values
- Returns:
Self
- Return type:
self
- sp_qc_engine(ids=None, program='xtb', model=None, keywords=None, local_options=None)[source]#
Calculate conformer energies with QCEngine interface.
- Parameters:
ids (Sequence[int] | None) – Conformer indices (1-indexed). If None, all are calculated.
program (str) – QCEngine program
model (dict[str, Any] | None) – QCEngine model
keywords (dict[str, Any] | None) – QCEngine keywords
local_options (dict[str, Any] | None) – QCEngine local options
- Returns:
Self
- Return type:
self
- Raises:
Exception – When trying to use RDKit with formal charges.
- write_xyz(file, ids=None, unit='kcal/mol', relative=True, separate=False, n_decimals=3)[source]#
Write conformers to xyz file.
- Parameters:
file (str | PathLike) – Filename or path object. Needs filename if separate=True
ids (Iterable[int] | None) – Conformer indices (1-indexed)
unit (str) – Output unit for energies in xyz comment field: ‘hartree’, ‘kcal/mol’, ‘kJ/mol’
relative (bool) – Whether to give energies relative to lowest energy conformer
separate (bool) – Whether to write conformers to separate xyz files
n_decimals (int) – Number of decimals for energies
- Raises:
TypeError – When separate=True and file is not str
- Return type:
None
- morfeus.conformer.boltzmann_average_dT(properties, energies, temperature=298.15)[source]#
Return the derivative of the Boltzmann average.
- Parameters:
properties (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – Conformer properties
energies (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – Conformer energies (a.u.)
temperature (float) – Temperature (K)
- Returns:
Derivative of Boltzmann average.
- Return type:
derivative
- morfeus.conformer.boltzmann_statistic(properties, energies, temperature=298.15, statistic='avg')[source]#
Compute Boltzmann statistic.
- Parameters:
properties (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – Conformer properties
energies (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – Conformer energies (a.u.)
temperature (float) – Temperature (K)
statistic (str) – Statistic to compute: ‘avg’, ‘var’ or ‘std’
- Returns:
Boltzmann statistic
- Return type:
result
- Raises:
ValueError – When statistic not specified correctly
- morfeus.conformer.boltzmann_weights(energies, temperature=298.15)[source]#
Compute Boltzmann weights.
- Parameters:
energies (_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]) – Conformer energies (a.u.)
temperature (float) – Temperature (K)
- Returns:
Conformer weights (normalized to unity)
- Return type:
weights
- morfeus.conformer.cli(smiles, generator='rdkit')[source]#
CLI for cone angle.
- Parameters:
smiles (str) – SMILES string
generator (str) – Confomer generator: ‘ob-ff’, ‘ob-ga’ or ‘rdkit’
- Returns:
Partially instantiated class
- Return type:
Any
- morfeus.conformer.conformers_from_ob_ff(mol, num_conformers=30, ff='MMFF94', method='systematic', rings=False)[source]#
Generates conformers based on the force field algorithm in OpenBabel.
Follows the recipe of the command line script obabel –conformer: https://github.com/openbabel/openbabel/blob/master/src/ops/conformer.cpp If an OBMol object with 3D coordinates is given, the conformer search will start from that structure.
- Parameters:
mol (str | ob.OBMol) – Molecule either as SMILES string or OBMol object.
num_conformers (int) – Maximum number of conformers
ff (str) – Force field supported by OpenBabel
method (str) – ‘fast’, random’, ‘systematic’ or ‘weighted’
rings (bool) – Sample ring torsions.
- Returns:
Elements as atomic numbers conformer_coordinates: Conformer coordinates (Å) connectivity_matrix: Connectivity matrix charges: Formal charges mol: OBMol object
- Return type:
elements
- morfeus.conformer.conformers_from_ob_ga(mol, num_conformers=None, num_children=None, mutability=None, convergence=None, score='rmsd', filter_method='steric', cutoff=0.8, vdw_factor=0.5, check_hydrogens=True)[source]#
Generates conformers based on the genetic algorithm in OpenBabel.
Follows the recipe of the command line script obabel –conformer: https://github.com/openbabel/openbabel/blob/master/src/ops/conformer.cpp If an OBMol object with 3D coordinates is given, the conformer search will start from that structure.
- Parameters:
mol (str | ob.OBMol) – Molecule either as SMILES string or OBMol object
num_conformers (int | None) – Maximum number of conformers
num_children (int | None) – Number of children to generate for each parent
mutability (float | None) – Mutation frequency
convergence (int | None) – Number of identical generations before convergence is reached
score (str) – Scoring function: ‘rmsd’, ‘min_rmsd’, ‘energy’, ‘min_energy’
filter_method (str) – Filtering algorithm: ‘steric’
cutoff (float) – Absolute distance in Ånström below which atoms are considered to clash
vdw_factor (float) – Scale factor applied to van der Waals radii for detecting clashes
check_hydrogens (bool) – Detect clashes with hydrogen atoms
- Returns:
Elements as atomic numbers conformer_coordinates: Conformer coordinates (Å) connectivity_matrix: Connectivity matrix charges: Formal charges mol: OBMol object
- Return type:
elements
- morfeus.conformer.conformers_from_rdkit(mol, n_confs=None, optimize='MMFF94', version=2, small_rings=True, macrocycles=True, random_seed=None, rmsd_thres=0.35, rmsd_symmetry=False, n_threads=1)[source]#
Generates conformers for an RDKit mol object.
Recipe based on J. Chem. Inf. Modeling 2012, 52, 1146.
- Parameters:
mol (str | Chem.Mol) – Molecule either as SMILES string or RDKit Mol object.
n_confs (int | None) – Number of conformers to generate. If None, a reasonable number will be set depending on the number of rotatable bonds.
optimize (str | None) – Force field used for conformer optimization: ‘MMFF94’, ‘MMFF94s’ or ‘UFF’. If None, conformers are not optimized.
version (int) – Version of the experimental torsion-angle preferences
small_rings (bool) – Whether to impose small ring torsion angle preferences
macrocycles (bool) – Whether to mpose macrocycle torsion angle preferences
random_seed (int | None) – Random seed for conformer generation
rmsd_thres (float | None) – Pruning RMSD threshold (Å)
rmsd_symmetry (bool) – Whether to use symmetry for RMSD pruning
n_threads (int) – Number of threads
- Returns:
Atomic symbols conformer_coordinates: Coordinates for all conformers (Å) energies: Conformer energies (a.u.) connectivity_matrix: Connectivity matrix with bond orders charges: Formal charges mol: RDKit Mol object. Only returned when return_mol=True
- Return type:
elements
- Raises:
Exception – When force field not found