arc.species.conformers

A module for (non-TS) species conformer generation

Note

variables that contain atom indices such as torsions and tops are 1-indexed, while atoms in Molecules are 0-indexed.

conformers is a list of dictionaries, each with the following keys:

{'xyz': <dict>,
 'index': <int>,
 'FF energy': <float>,
 'source': <str>,
 'torsion_dihedrals': {<torsion tuple 0>: angle 0,
                       <torsion tuple 1>: angle 1,
 }

Module workflow:

generate_conformers
    generate_force_field_conformers
        get_force_field_energies, rdkit_force_field or openbabel_force_field_on_rdkit_conformers,
        determine_dihedrals
    deduce_new_conformers
        get_torsion_angles, determine_torsion_symmetry, determine_torsion_sampling_points,
        change_dihedrals_and_force_field_it
    get_lowest_confs
arc.species.conformers.change_dihedrals_and_force_field_it(label, mol, xyz, torsions, new_dihedrals, optimize=True, force_field='MMFF94s')[source]

Change dihedrals of specified torsions according to the new dihedrals specified, and get FF energies.

Example:

torsions = [(1, 2, 3, 4), (9, 4, 7, 1)]
new_dihedrals = [[90, 120], [90, 300], [180, 270], [30, 270]]

This will calculate the energy of the original conformer (defined using xyz). We iterate through new_dihedrals. The torsions are set accordingly and the energy and xyz of the newly generated conformer are kept. We assume that each list entry in new_dihedrals is of the length of the torsions list (2 in the example).

Parameters:
  • label (str) – The species’ label.

  • mol (Molecule) – The RMG molecule with the connectivity information.

  • xyz (dict) – The base 3D geometry to be changed.

  • torsions (list) – Entries are torsion tuples for which the dihedral will be changed relative to xyz.

  • new_dihedrals (list) – Entries are same size lists of dihedral angles (floats) corresponding to the torsions.

  • optimize (bool, optional) – Whether to optimize the coordinates using FF. True to optimize.

  • force_field (str, optional) – The type of force field to use.

Returns:

  • The conformer FF energies corresponding to the list of dihedrals.

  • The conformer xyz geometries corresponding to the list of dihedrals.

Return type:

Tuple[list, list]

arc.species.conformers.cheat_sheet(mol_list: Union[List[Molecule], Molecule]) Optional[List[Dict]][source]

Check if the species is in the cheat sheet, and return its correct xyz if it is.

Parameters:

mol_list (Union[List[Molecule], Molecule]) – Molecule objects to consider (or Molecule, resonance structures will be generated).

Returns: Optional[lis[Dict]]

arc.species.conformers.check_special_non_rotor_cases(mol, top1, top2)[source]

Check whether one of the tops correspond to a special case which does not have a torsional mode. Checking for R-[C,N]#[N,[CH],[C]] groups, such as: in cyano groups (R-C#N`), C#C groups (R-C#CH or R-C#[C]), and azide groups: (R-N#N).

Parameters:
  • mol (Molecule) – The RMG molecule.

  • top1 (list) – Entries are atom indices (1-indexed) on one side of the torsion, inc. one of the pivotal atoms.

  • top2 (list) – Entries are atom indices (1-indexed) on the other side of the torsion, inc. the other pivotal atom.

Returns:

True if this is indeed a special case which should not be treated as a torsional mode.

Return type:

bool

arc.species.conformers.chirality_dict_to_tuple(chirality_dict)[source]

A helper function for using the chirality dictionary of a conformer as a key in the enantiomers_dict by converting it to a tuple deterministically.

Parameters:

chirality_dict (dict) – The chirality dictionary of a conformer.

Raises:

ConformerError – If the chirality values are wrong.

Returns:

A deterministic tuple representation of the chirality dictionary.

Return type:

tuple

arc.species.conformers.conformers_combinations_by_lowest_conformer(label, mol, base_xyz, multiple_tors, multiple_sampling_points, len_conformers=- 1, force_field='MMFF94s', max_combination_iterations=25, torsion_angles=None, multiple_sampling_points_dict=None, wells_dict=None, de_threshold=None, plot_path=False, symmetries=None)[source]

Iteratively modify dihedrals in the lowest conformer (each iteration deduces a new lowest conformer), until convergence.

Parameters:
  • label (str) – The species’ label.

  • mol (Molecule) – The RMG molecule with the connectivity information.

  • base_xyz (dict) – The base 3D geometry to be changed.

  • multiple_tors (list) – Entries are torsion tuples of non-symmetric torsions.

  • multiple_sampling_points (list) – Entries are lists of dihedral angles (sampling points), respectively correspond to torsions in multiple_tors.

  • len_conformers (int, optional) – The length of the existing conformers list (for consecutive numbering).

  • de_threshold (float, optional) – An energy threshold (in kJ/mol) above which wells in a torsion will not be considered.

  • force_field (str, optional) – The type of force field to use.

  • max_combination_iterations (int, optional) – The max num of times to iteratively search for the lowest conformer.

  • torsion_angles (dict, optional) – The torsion angles. Keys are torsion tuples, values are lists of all corresponding angles from conformers.

  • multiple_sampling_points_dict (dict, optional) – Keys are torsion tuples, values are respective sampling points.

  • wells_dict (dict, optional) – Keys are torsion tuples, values are well dictionaries.

  • plot_path (str, optional) – A folder path in which the plot will be saved. If None, the plot will not be shown (nor saved).

  • symmetries (dict, optional) – Keys are tuples scan indices (1-indexed), values are internal rotation symmetry numbers (sigma).

Returns:

New conformer combinations, entries are conformer dictionaries.

Return type:

list

arc.species.conformers.deduce_new_conformers(label, conformers, torsions, tops, mol_list, smeared_scan_res=None, plot_path=None, combination_threshold=1000, force_field='MMFF94s', max_combination_iterations=25, diastereomers=None, de_threshold=None)[source]

By knowing the existing torsion wells, get the geometries of all important conformers. Validate that atoms don’t collide in the generated conformers (don’t consider ones where they do).

Parameters:
  • label (str) – The species’ label.

  • conformers (list) – Entries are conformer dictionaries.

  • torsions (list) – A list of all possible torsion angles in the molecule, each torsion angles list is sorted.

  • tops (list) – A list of tops corresponding to torsions.

  • mol_list (list) – A list of RMG Molecule objects.

  • smeared_scan_res (float, optional) – The resolution (in degrees) for scanning smeared wells.

  • plot_path (str, optional) – A folder path in which the plot will be saved. If None, the plot will not be shown (nor saved).

  • combination_threshold (int, optional) – A threshold below which all combinations will be generated.

  • force_field (str, optional) – The type of force field to use.

  • max_combination_iterations (int, optional) – The max num of times to iteratively search for the lowest conformer.

  • diastereomers (list, optional) – Entries are xyz’s in a dictionary format or conformer structures representing specific diastereomers to keep.

  • de_threshold (float, optional) – An energy threshold (in kJ/mol) above which wells in a torsion will not be considered.

Returns:

  • The deduced conformers.

  • Keys are torsion tuples.

Return type:

Tuple[list, dict]

arc.species.conformers.determine_chirality(conformers, label, mol, force=False)[source]

Determines the Cahn–Ingold–Prelog (CIP) chirality (R or S) of atoms in the conformer, as well as the CIP chirality of double bonds (E or Z).

Parameters:
  • conformers (list) – Entries are conformer dictionaries.

  • label (str) – The species’ label.

  • mol (RMG Molecule or RDKit RDMol) – The molecule object with connectivity and bond order information.

  • force (bool, optional) – Whether to override data, True to override, default is False.

Returns:

Conformer dictionaries with updated with ‘chirality’. conformer['chirality'] is a dictionary.

Keys are either a 1-length tuple of atom indices (for chiral atom centers) or a 2-length tuple of atom indices (for chiral double bonds), values are either ‘R’ or ‘S’ for chiral atom centers (or ‘NR’ or ‘NS’ for chiral nitrogen centers), or ‘E’ or ‘Z’ for chiral double bonds. All atom indices are 0-indexed.

Return type:

list

arc.species.conformers.determine_dihedrals(conformers, torsions)[source]

For each conformer in conformers determine the respective dihedrals.

Parameters:
  • conformers (list) – Entries are conformer dictionaries.

  • torsions (list) – All possible torsions in the molecule.

Returns:

Entries are conformer dictionaries.

Return type:

list

arc.species.conformers.determine_number_of_conformers_to_generate(label: str, heavy_atoms: int, torsion_num: int, mol: Optional[Molecule] = None, xyz: Optional[dict] = None, minimalist: bool = False) Tuple[int, int][source]

Determine the number of conformers to generate using molecular mechanics.

Parameters:
  • label (str) – The species’ label.

  • heavy_atoms (int) – The number of heavy atoms in the molecule.

  • torsion_num (int) – The number of potential torsions in the molecule.

  • mol (Molecule, optional) – The RMG Molecule object.

  • xyz (dict, optional) – The xyz coordinates.

  • minimalist (bool, optional) – Whether to return a small number of conformers, useful when this is just a guess before fitting a force field. True to be minimalistic.

Raises:

ConformerError – If the number of conformers to generate cannot be determined.

Returns:

  • The number of conformers to generate.

  • The number of chiral centers.

Return type:

Tuple[int, int]

arc.species.conformers.determine_rotors(mol_list)[source]

Determine possible unique rotors in the species to be treated as hindered rotors.

Parameters:

mol_list (list) – Localized structures (Molecule objects) by which all rotors will be determined.

Returns:

  • A list of indices of scan pivots.

  • A list of indices of top atoms (including one of the pivotal atoms) corresponding to the torsions.

Return type:

Tuple[list, list]

arc.species.conformers.determine_smallest_atom_index_in_scan(atom1: Atom, atom2: Atom, mol: Molecule) int[source]

Determine the smallest atom index in mol connected to atom1 which is not atom2. Returns a heavy atom if available, otherwise a hydrogen atom. Useful for deterministically determining the indices of four atom in a scan. This function assumes there ARE additional atoms connected to atom1, and that atom2 is not a hydrogen atom.

Parameters:
  • atom1 (Atom) – The atom whose neighbors will be searched.

  • atom2 (Atom) – An atom connected to atom1 to exclude (a pivotal atom).

  • mol (Molecule) – The molecule to process.

Returns:

The smallest atom index (1-indexed) connected to atom1 which is not atom2.

Return type:

int

arc.species.conformers.determine_torsion_sampling_points(label, torsion_angles, smeared_scan_res=None, symmetry=1)[source]

Determine how many points to consider in each well of a torsion for conformer combinations.

Parameters:
  • label (str) – The species’ label.

  • torsion_angles (list) – Well angles in the torsion.

  • smeared_scan_res (float, optional) – The resolution (in degrees) for scanning smeared wells.

  • symmetry (int, optional) – The torsion symmetry number.

Returns:

  • Sampling points for the torsion.

  • Each entry is a well dictionary with the keys start_idx, end_idx, start_angle, end_angle, angles.

Return type:

Tuple[list, list]

arc.species.conformers.determine_torsion_symmetry(label, top1, mol_list, torsion_scan)[source]

Check whether a torsion is symmetric. If a torsion well is “well defined” and not smeared, it could be symmetric. Check the groups attached to the rotor pivots to determine whether it is indeed symmetric We don’t care about the actual rotor symmetry number here, since we plan to just use the first well (they’re all the same).

Parameters:
  • label (str) – The species’ label.

  • top1 (list) – A list of atom indices on one side of the torsion, including the pivotal atom.

  • mol_list (list) – A list of molecules.

  • torsion_scan (list) – The angles corresponding to this torsion from all conformers.

Returns:

The rotor symmetry number.

Return type:

int

arc.species.conformers.determine_well_width_tolerance(mean_width)[source]

Determine the tolerance by which well widths are determined to be nearly equal.

Fitted to a polynomial trend line for the following data of (mean, tolerance) pairs:

(100, 0.11), (60, 0.13), (50, 0.15), (25, 0.25), (5, 0.50), (1, 0.59)
Parameters:

mean_width (float) – The mean well width in degrees.

Returns:

The tolerance.

Return type:

float

arc.species.conformers.embed_rdkit(label, mol, num_confs=None, xyz=None)[source]

Generate unoptimized conformers in RDKit. If xyz is not given, random conformers will be generated.

Parameters:
  • label (str) – The species’ label.

  • mol (RMG Molecule or RDKit RDMol) – The molecule object with connectivity and bond order information.

  • num_confs (int, optional) – The number of random 3D conformations to generate.

  • xyz (dict, optional) – The 3D coordinates.

Returns:

An RDKIt molecule with embedded conformers.

Return type:

Optional[RDMol]

arc.species.conformers.find_internal_rotors(mol)[source]

Locates the sets of indices corresponding to every internal rotor (1-indexed).

Parameters:

mol (Molecule) – The molecule for which rotors will be determined.

Returns:

Entries are rotor dictionaries with the four-atom scan coordinates, the pivots, and the smallest top.

Return type:

list

arc.species.conformers.generate_all_combinations(label, mol, base_xyz, multiple_tors, multiple_sampling_points, len_conformers=- 1, torsions=None, force_field='MMFF94s')[source]

Generate all combinations of torsion wells from a base conformer.

Parameters:
  • label (str) – The species’ label.

  • mol (Molecule) – The RMG molecule with the connectivity information.

  • base_xyz (dict) – The base 3D geometry to be changed.

  • multiple_tors (list) – Entries are torsion tuples of non-symmetric torsions.

  • multiple_sampling_points (list) – Entries are lists of dihedral angles (sampling points), respectively correspond to torsions in multiple_tors.

  • len_conformers (int, optional) – The length of the existing conformers list (for consecutive numbering).

  • force_field (str, optional) – The type of force field to use.

  • torsions (list, optional) – A list of all possible torsions in the molecule. Will be determined if not given.

Returns:

New conformer combinations, entries are conformer dictionaries.

Return type:

list

arc.species.conformers.generate_conformer_combinations(label, mol, base_xyz, hypothetical_num_comb, multiple_tors, multiple_sampling_points, combination_threshold=1000, len_conformers=- 1, force_field='MMFF94s', max_combination_iterations=25, plot_path=None, torsion_angles=None, multiple_sampling_points_dict=None, wells_dict=None, de_threshold=None, symmetries=None)[source]

Call either conformers_combinations_by_lowest_conformer() or generate_all_combinations(), according to the hypothetical_num_comb.

Parameters:
  • label (str) – The species’ label.

  • mol (Molecule) – The RMG molecule with the connectivity information.

  • base_xyz (dict) – The base 3D geometry to be changed.

  • hypothetical_num_comb (int) – The number of combinations that could be generated by changing dihedrals, considering symmetry but not considering atom collisions.

  • combination_threshold (int, optional) – A threshold below which all combinations will be generated.

  • multiple_tors (list) – Entries are torsion tuples of non-symmetric torsions.

  • multiple_sampling_points (list) – Entries are lists of dihedral angles (sampling points), respectively correspond to torsions in multiple_tors.

  • len_conformers (int, optional) – The length of the existing conformers list (for consecutive numbering).

  • de_threshold (float, optional) – An energy threshold (in kJ/mol) above which wells in a torsion will not be considered.

  • force_field (str, optional) – The type of force field to use.

  • max_combination_iterations (int, optional) – The max num of times to iteratively search for the lowest conformer.

  • torsion_angles (dict, optional) – The torsion angles. Keys are torsion tuples, values are lists of all corresponding angles from conformers.

  • multiple_sampling_points_dict (dict, optional) – Keys are torsion tuples, values are respective sampling points.

  • wells_dict (dict, optional) – Keys are torsion tuples, values are well dictionaries.

  • plot_path (str, optional) – A folder path in which the plot will be saved. If None, the plot will not be shown (nor saved).

  • symmetries (dict, optional) – Keys are tuples scan indices (1-indexed), values are internal rotation symmetry numbers (sigma).

Returns:

New conformer combinations, entries are conformer dictionaries.

Return type:

list

arc.species.conformers.generate_conformers(mol_list: Union[List[Molecule], Molecule], label, xyzs=None, torsions=None, tops=None, charge=0, multiplicity=None, num_confs_to_generate=None, n_confs=10, e_confs=5.0, de_threshold=None, smeared_scan_res=None, combination_threshold=None, force_field='MMFF94s', max_combination_iterations=None, diastereomers=None, return_all_conformers=False, plot_path=None, print_logs=True) Optional[Union[list, Tuple[list, list]]][source]

Generate conformers for (non-TS) species starting from a list of RMG Molecules. (resonance structures are assumed to have already been generated and included in the molecule list)

Parameters:
  • mol_list (Union[List[Molecule], Molecule]) – Molecule objects to consider (or Molecule, resonance structures will be generated).

  • label (str) – The species’ label.

  • xyzs (list) – A list of user guess xyzs that will also be taken into account, each in a dict format.

  • torsions (list, optional) – A list of all possible torsions in the molecule. Will be determined if not given.

  • tops (list, optional) – A list of tops corresponding to torsions. Will be determined if not given.

  • charge (int, optional) – The species charge. Used to perceive a molecule from xyz.

  • multiplicity (int, optional) – The species multiplicity. Used to perceive a molecule from xyz.

  • num_confs_to_generate (int, optional) – The number of conformers to generate (can be determined automatically)

  • n_confs (int, optional) – The number of conformers to return.

  • e_confs (float, optional) – The energy threshold in kJ/mol above the lowest energy conformer below which all (unique) generated conformers will be returned.

  • de_threshold (float, optional) – Energy threshold (in kJ/mol) above which wells will not be considered.

  • smeared_scan_res (float, optional) – The resolution (in degrees) for scanning smeared wells.

  • combination_threshold (int, optional) – A threshold below which all combinations will be generated.

  • force_field (str, optional) – The type of force field to use (MMFF94, MMFF94s, UFF, GAFF).

  • max_combination_iterations (int, optional) – The maximum number of times to iteratively search for the lowest conformer.

  • diastereomers (list, optional) – Entries are xyz’s in a dictionary format or conformer structures representing specific diastereomers to keep.

  • return_all_conformers (bool, optional) – Whether to return the full conformers list of conformer dictionaries In addition to the lowest conformers list. Tru to return it.

  • plot_path (str, optional) – A folder path in which the plot will be saved. If None, the plot will not be shown (nor saved).

  • print_logs (bool, optional) – Whether define a logger so logs are also printed to stdout. Useful when run outside of ARC. True to print.

Raises:
  • ConformerError – If something goes wrong.

  • TypeError – If xyzs has entries of a wrong type.

Returns: Optional[Union[list, Tuple[list, list]]]
  • Lowest conformers

  • Lowest conformers and all new conformers.

arc.species.conformers.generate_diatomic_conformer(symbol_1: str, symbol_2: str, multiplicity: Optional[int] = None) dict[source]

Generate a conformer for a diatomic species. Data from CCCBDB.

Parameters:
  • symbol_1 (str) – The atomic symbol of atom 1.

  • symbol_2 (str) – The atomic symbol of atom 2.

  • multiplicity (int, optional) – The diatomic species multiplicity

Returns:

The diatomic conformer.

Return type:

dict

arc.species.conformers.generate_force_field_conformers(label, mol_list, torsion_num, charge, multiplicity, xyzs=None, num_confs=None, force_field='MMFF94s') List[dict][source]

Generate conformers using RDKit and OpenBabel and optimize them using a force field Also consider user guesses in xyzs.

Parameters:
  • label (str) – The species’ label.

  • mol_list (list) – Entries are Molecule objects representing resonance structures of a chemical species.

  • xyzs (list, optional) – Entries are xyz coordinates in dict format, given as initial guesses.

  • torsion_num (int) – The number of torsions identified in the molecule.

  • charge (int) – The net charge of the species.

  • multiplicity (int) – The species spin multiplicity.

  • num_confs (int, optional) – The number of conformers to generate.

  • force_field (str, optional) – The type of force field to use.

Raises:

ConformerError – If xyzs is given and it is not a list, or its entries are not strings.

Returns:

Entries are conformer dictionaries.

Return type:

list

arc.species.conformers.generate_monoatomic_conformer(symbol: str) dict[source]

Generate a conformer for a monoatomic species.

Parameters:

symbol (str) – The atomic symbol.

Returns:

The monoatomic conformer.

Return type:

dict

arc.species.conformers.get_force_field_energies(label: str, mol: Molecule, num_confs: Optional[int] = None, xyz: Optional[dict] = None, force_field: str = 'MMFF94s', try_uff: bool = True, optimize: bool = True, try_ob: bool = False, suppress_warning: bool = False) Tuple[list, list][source]

Determine force field energies using RDKit. If num_confs is given, random 3D geometries will be generated. If xyz is given, it will be directly used instead. The coordinates are returned in the order of atoms in mol.

Parameters:
  • label (str) – The species’ label.

  • mol (Molecule) – The RMG molecule object with connectivity and bond order information.

  • num_confs (int, optional) – The number of random 3D conformations to generate.

  • xyz (dict, optional) – The 3D coordinates guess.

  • force_field (str, optional) – The type of force field to use.

  • try_uff (bool, optional) – Whether to try UFF if MMFF94(s) fails. True by default.

  • optimize (bool, optional) – Whether to first optimize the conformer using FF. True to optimize.

  • try_ob (bool, optional) – Whether to try OpenBabel if RDKit fails. True to try, True by default.

  • suppress_warning (bool, optional) – Whether to suppress OpenBabel warnings. False by default.

Raises:

ConformerError – If conformers could not be generated.

Returns:

  • Entries are xyz coordinates, each in a dict format.

  • Entries are the FF energies (in kJ/mol).

Return type:

Tuple[list, list]

arc.species.conformers.get_lowest_confs(label: str, confs: Union[dict, list], n: Optional[int] = 10, e: Optional[float] = 5.0, energy: str = 'FF energy') list[source]

Get the most stable conformer.

Parameters:
  • label (str) – The species’ label.

  • confs (dict, list) – Entries are either conformer dictionaries or a length two list of xyz coordinates and energy

  • n (int, optional) – Number of lowest conformers to return.

  • e (float, optional) – The energy threshold above the lowest energy conformer in kJ/mol below which all conformers will be returned.

  • energy (str, optional) – The energy attribute to search by. Currently only ‘FF energy’ is supported.

Raises:

ConformerError – If n < 1, e < 0, both n and e are None, or if no conformers are given.

Returns:

Conformer dictionaries.

Return type:

list

arc.species.conformers.get_lowest_diastereomers(label, mol, conformers, diastereomers=None)[source]

Get the 2^(n-1) diastereomers with the lowest energy (where n is the number of chiral centers in the molecule). We exclude enantiomers (mirror images where ALL chiral centers simultaneously invert). If a specific diastereomer is given (in an xyz dict form), then only the lowest conformer with the same chirality will be returned.

Parameters:
  • label (str) – The species’ label.

  • mol (Molecule) – The 2D graph representation of the molecule.

  • conformers (list) – Entries are conformer dictionaries.

  • diastereomers (list, optional) – Entries are xyz’s in a dictionary format or conformer structures representing specific diastereomers to keep.

Raises:

ConformerError – If diastereomers is not None and is of wrong type, or if conformers with the requested chirality combination could not be generated.

Returns:

Entries are lowest energy diastereomeric conformer dictionaries to consider.

Return type:

list

arc.species.conformers.get_number_of_chiral_centers(label, mol, conformer=None, xyz=None, just_get_the_number=True)[source]

Determine the number of chiral centers by type. Either conformer or xyz must be given.

Parameters:
  • label (str) – The species label.

  • mol (Molecule) – The RMG Molecule object.

  • conformer (dict, optional) – A conformer dictionary.

  • xyz (dict, optional) – The xyz coordinates.

  • just_get_the_number (bool, optional) – Return the number of chiral centers regardless of their type.

Raises:

InputError – If neither conformer nor xyz were given.

Returns:

Keys are types of chiral sites (‘C’ for carbon, ‘N’ for nitrogen, ‘D’ for double bond), values are the number of chiral centers of each type. If just_get_the_number is True, just returns the number of chiral centers (integer).

Return type:

Optional[dict, int]

arc.species.conformers.get_top_element_count(mol, top)[source]

Returns the element count for the molecule considering only the atom indices in top.

Parameters:
  • mol (Molecule) – The molecule to consider.

  • top (list) – The atom indices to consider.

Returns:

The element count, keys are tuples of (element symbol, isotope number), values are counts.

Return type:

dict

arc.species.conformers.get_torsion_angles(label, conformers, torsions)[source]

Populate each torsion pivots with all available angles from the generated conformers.

Parameters:
  • label (str) – The species’ label.

  • conformers (list) – The conformers from which to extract the angles.

  • torsions (list) – The torsions to consider.

Returns:

The torsion angles. Keys are torsion tuples, values are lists of all corresponding angles from conformers.

Return type:

dict

arc.species.conformers.get_wells(label, angles, blank=20)[source]

Determine the distinct wells from a list of angles.

Parameters:
  • label (str) – The species’ label.

  • angles (list) – The angles in the torsion.

  • blank (int, optional) – The blank space between wells.

Returns:

Entry are well dicts with keys: start_idx, end_idx, start_angle, end_angle, angles.

Return type:

list

arc.species.conformers.identify_chiral_nitrogen_centers(mol)[source]

Identify the atom indices corresponding to a chiral nitrogen centers in a molecule (umbrella modes).

Parameters:

mol (Molecule) – The molecule to be analyzed.

Raises:

TypeError – If mol is of wrong type.

Returns:

Atom numbers (0-indexed) representing chiral nitrogen centers in the molecule (umbrella modes).

Return type:

list

arc.species.conformers.initialize_log(verbose=20)[source]

Set up a simple logger for stdout printing (not saving into as log file).

Parameters:

verbose (int, optional) – Specify the amount of log text seen.

arc.species.conformers.inverse_chirality_symbol(symbol)[source]

Inverses a chirality symbol, e.g., the ‘R’ character to ‘S’, or ‘NS’ to ‘NR’. Note that chiral double bonds (‘E’ and ‘Z’) must not be inversed (they are not mirror images of each other).

Parameters:

symbol (str) – The chirality symbol.

Raises:

InputError – If symbol could not be recognized.

Returns:

The inverse chirality symbol.

Return type:

str

arc.species.conformers.mix_rdkit_and_openbabel_force_field(label, mol, num_confs=None, xyz=None, force_field='GAFF', try_ob=False)[source]

Optimize conformers using a force field (GAFF, MMFF94s, MMFF94, UFF, Ghemical) Use RDKit to generate the random conformers (OpenBabel isn’t good enough), but use OpenBabel to optimize them (RDKit doesn’t have GAFF).

Parameters:
  • label (str) – The species’ label.

  • mol (Molecule, optional) – The RMG molecule object with connectivity and bond order information.

  • num_confs (int, optional) – The number of random 3D conformations to generate.

  • xyz (string or list, optional) – The 3D coordinates in either a string or an array format.

  • force_field (str, optional) – The type of force field to use.

  • try_ob (bool, optional) – Whether to try OpenBabel if RDKit fails. True to try, False by default.

Returns:

  • Entries are optimized xyz’s in a list format.

  • Entries are float numbers representing the energies in kJ/mol.

Return type:

Tuple[list, list]

arc.species.conformers.openbabel_force_field(label, mol, num_confs=None, xyz=None, force_field='GAFF', method='diverse')[source]

Optimize conformers using a force field (GAFF, MMFF94s, MMFF94, UFF, Ghemical).

Parameters:
  • label (str) – The species’ label.

  • mol (Molecule, optional) – The RMG molecule object with connectivity and bond order information.

  • num_confs (int, optional) – The number of random 3D conformations to generate.

  • xyz (dict, optional) – The 3D coordinates.

  • force_field (str, optional) – The type of force field to use.

  • method (str, optional) – The conformer searching method to use in OpenBabel. For method description, see https://openbabel.org/dev-api/group__conformer.shtml

Returns:

  • Entries are optimized xyz’s in a list format.

  • Entries are float numbers representing the energies in kJ/mol.

Return type:

Tuple[list, list]

arc.species.conformers.openbabel_force_field_on_rdkit_conformers(label, rd_mol, force_field='MMFF94s', optimize=True)[source]

Optimize RDKit conformers by OpenBabel using a force field (MMFF94 or MMFF94s are recommended). This is a fallback method when RDKit fails to generate force field optimized conformers.

Parameters:
  • label (str) – The species’ label.

  • rd_mol (RDKit RDMol) – The RDKit molecule with embedded conformers to optimize.

  • force_field (str, optional) – The type of force field to use.

  • optimize (bool, optional) – Whether to first optimize the conformer using FF. True to optimize.

Returns:

  • Entries are optimized xyz’s in a dictionary format.

  • Entries are float numbers representing the energies (in kJ/mol).

Return type:

Tuple[list, list]

arc.species.conformers.prune_enantiomers_dict(label, enantiomers_dict)[source]

A helper function for screening out enantiomers from the enantiomers_dict, leaving only diastereomers (so removing all exact mirror images). Note that double bond chiralities ‘E’ and ‘Z’ are not mirror images of each other, and are not pruned out.

Parameters:
  • label (str) – The species’ label.

  • enantiomers_dict (dict) – Keys are chirality tuples, values are conformer structures.

Returns:

The pruned enantiomers_dict.

Return type:

dict

arc.species.conformers.rdkit_force_field(label: str, rd_mol: Optional[EditableMol], mol: Optional[Molecule] = None, num_confs: Optional[int] = None, force_field: str = 'MMFF94s', try_uff: bool = True, optimize: bool = True, try_ob: bool = False) Tuple[list, list][source]

Optimize RDKit conformers using a force field (MMFF94 or MMFF94s are recommended). For UFF see: https://www.rdkit.org/docs/source/rdkit.Chem.rdForceFieldHelpers.html#rdkit.Chem.rdForceFieldHelpers.UFFOptimizeMoleculeConfs

Parameters:
  • label (str) – The species’ label.

  • rd_mol (RDKit RDMol) – The RDKit molecule with embedded conformers to optimize.

  • mol (Molecule, optional) – The RMG molecule object with connectivity and bond order information.

  • num_confs (int, optional) – The number of random 3D conformations to generate.

  • force_field (str, optional) – The type of force field to use.

  • try_uff (bool, optional) – Whether to try UFF if MMFF94(s) fails. True by default.

  • optimize (bool, optional) – Whether to first optimize the conformer using FF. True to optimize.

  • try_ob (bool, optional) – Whether to try OpenBabel if RDKit fails. True to try, False by default.

Returns:

  • Entries are optimized xyz’s in a dictionary format.

  • Entries are float numbers representing the energies.

Return type:

Tuple[list, list]

arc.species.conformers.read_rdkit_embedded_conformer_i(rd_mol, i, rd_index_map=None)[source]

Read coordinates from RDKit conformers.

Parameters:
  • rd_mol (RDKit RDMol) – The RDKit molecule with embedded conformers to optimize.

  • i (int) – The conformer index from rd_mol to read.

  • rd_index_map (list, optional) – An atom map dictionary to reorder the xyz. Keys are rdkit atom indices, values are RMG mol atom indices.

Returns:

xyz coordinates.

Return type:

dict

arc.species.conformers.read_rdkit_embedded_conformers(label, rd_mol, i=None, rd_index_map=None)[source]

Read coordinates from RDKit conformers.

Parameters:
  • label (str) – The species’ label.

  • rd_mol (RDKit RDMol) – The RDKit molecule with embedded conformers to optimize.

  • i (int, optional) – The conformer index from rd_mol to read. If None, all will be read,

  • rd_index_map (list, optional) – An atom map dictionary to reorder the xyz. Requires mol to not be None.

Returns:

Entries are xyz coordinate dicts.

Return type:

list

arc.species.conformers.replace_n_with_c_in_mol(mol, chiral_nitrogen_centers)[source]

Replace nitrogen atoms (pre-identified as chiral centers) with carbon atoms, replacing the lone electron pair (assuming just one exists) with a hydrogen or a halogen atom, preserving any radical electrons on the nitrogen atom.

Parameters:
  • mol (Molecule) – The molecule to be analyzed.

  • chiral_nitrogen_centers (list) – The 0-index of chiral (umbrella mode) nitrogen atoms in the molecule.

Raises:

ConformerError – If any of the atoms indicated by chiral_nitrogen_centers could not be a chiral nitrogen atom.

Returns:

  • A copy of the molecule with replaced N atoms.

  • Elements inserted in addition to the C atom, ordered as in chiral_nitrogen_centers.

Return type:

Tuple[Molecule, list]

arc.species.conformers.replace_n_with_c_in_xyz(label, mol, xyz, chiral_nitrogen_centers, elements_to_insert)[source]

Replace nitrogen atoms (pre-identified as chiral centers) with carbon atoms, replacing the lone electron pair (assuming just one exists) with a hydrogen or a halogen atom.

Parameters:
  • label (str) – The species label.

  • mol (Molecule) – The respective molecule object.

  • xyz (dict) – The 3D coordinates to process.

  • chiral_nitrogen_centers (list) – The 0-index of chiral (umbrella mode) nitrogen atoms in the molecule.

  • elements_to_insert (list) – The element (H/F/Cl/I) to insert in addition to C per nitrogen center.

Returns:

The coordinates with replaced N atoms.

Return type:

dict

arc.species.conformers.to_group(mol, atom_indices)[source]

This method converts a defined part of a Molecule into a Group.

Parameters:
  • mol (Molecule) – The base molecule.

  • atom_indices (list) – 0-indexed atom indices corresponding to atoms in mol to be included in the group.

Returns:

A group consisting of the desired atoms in mol.

Return type:

Group

arc.species.conformers.translate_group(mol, xyz, pivot, anchor, vector)[source]

Translate a group (a set of atoms from the pivot towards the anchor and onwards) by changing its pivot -> anchor vector to the desired new vector. Keep the relative distances between the group’s atoms constant, as well as the distance between the anchor and the vector atoms.

Parameters:
  • mol (Molecule) – The 2D graph representation of the molecule.

  • xyz (dict) – The 3D coordinates of the molecule with the same atom order as in mol.

  • pivot (int) – The 0-index of the pivotal atom around which groups are to be translated.

  • anchor (int) – The 0-index of an anchor atom. The group is defined from the pivot atom to the anchor atom, including all other atoms in the molecule connected to the anchor. The pivot and anchor atoms should not have another path connecting them such as a ring.

  • vector (list) – The new vector by which the group will be translated.

Returns:

The translated coordinates.

Return type:

dict

arc.species.conformers.translate_groups(label, mol, xyz, pivot)[source]

Exchange between two groups in a molecule. The groups cannot share a ring with the pivotal atom. The function does not change the atom order, just the coordinates of atoms. If the pivotal atom has exactly one lone pair, consider it as well as a dummy atom in translations.

Parameters:
  • label (str) – The species’ label.

  • mol (Molecule) – The 2D graph representation of the molecule.

  • xyz (dict) – A string-format 3d coordinates of the molecule with the same atom order as in mol.

  • pivot (int) – The 0-index of the pivotal atom around which groups are to be translated.

Returns:

The translated coordinates.

Return type:

dict

arc.species.conformers.update_mol(mol)[source]

Update atom types, multiplicity, and atom charges in the molecule.

Parameters:

mol (Molecule) – The molecule to update.

Returns:

the updated molecule.

Return type:

Molecule