arc.species.converter

A module for performing various species-related format conversions.

arc.species.converter.add_lone_pairs_by_atom_valance(mol)[source]

A helper function for assigning lone pairs instead of carbenes/nitrenes if not identified automatically, and they are missing according to the given multiplicity.

Parameters:

mol (Molecule) – The Molecule object to process.

arc.species.converter.add_rads_by_atom_valance(mol)[source]

A helper function for assigning radicals if not identified automatically, and they are missing according to the given multiplicity. We assume here that all partial charges were already set, but this assumption could be wrong. Note: This implementation might also be problematic for aromatic species with undefined bond orders.

Parameters:

mol (Molecule) – The Molecule object to process.

arc.species.converter.check_isomorphism(mol1, mol2, filter_structures=True, convert_to_single_bonds=False)[source]

Convert mol1 and mol2 to RMG Species objects, and generate resonance structures. Then check Species isomorphism. This function first makes copies of the molecules, since isIsomorphic() changes atom orders.

Parameters:
  • mol1 (Molecule) – An RMG Molecule object.

  • mol2 (Molecule) – An RMG Molecule object.

  • filter_structures (bool, optional) – Whether to apply the filtration algorithm when generating resonance structures. True to apply.

  • convert_to_single_bonds (bool, optional) – Whether to convert both molecules to single bonds, avoiding a bond order comparison (only compares connectivity). Resonance structures will not be generated.

Returns:

Whether one of the molecules in the Species derived from mol1

is isomorphic to one of the molecules in the Species derived from mol2. True if it is.

Return type:

bool

arc.species.converter.check_xyz_dict(xyz: Union[dict, str], project_directory: Optional[str] = None) Optional[dict][source]

Check that the xyz dictionary entered is valid. If it is a string, convert it. If it is a Z matrix, convert it to cartesian coordinates, If isotopes are not in xyz_dict, common values will be added. If a part of the xyz structure is a np.ndarray type, convert it by always calling xyz_from_data().

Parameters:
  • xyz (Union[dict, str]) – The xyz dictionary.

  • project_directory (str, optional) – The path to the project directory.

Raises:

ConverterError – If xyz is of wrong type or is missing symbols or coords.

Returns: Optional[dict]

The cartesian coordinates in a dictionary format.

arc.species.converter.check_zmat_dict(zmat: Union[dict, str]) dict[source]

Check that the zmat dictionary entered is valid. If it is a string, convert it. If it represents cartesian coordinates, convert it to internal coordinates. If a map isn’t given, a trivial one will be added.

Parameters:

zmat (dict, str) – The zmat dictionary.

Raises:

ConverterError – If zmat is of wrong type or is missing vars or coords.

arc.species.converter.cluster_confs_by_rmsd(xyzs: Iterable[Dict[str, tuple]], rmsd_threshold: float = 0.01) Tuple[Dict[str, tuple]][source]

Cluster conformers with the same atom orders using RMSD of distance matrices. Works for both TS and non-TS conformers. Intended for finding structurally distinct conformers from a pool of conformers. Suitable scenario: Filter a pool of conformers with their geometry optimized at some level. Not suitable for clustering conformers (not optimized) that are sampling of a well or a saddle point (these conformers may have large difference in RMSE, but they really should be representing the same well or saddle point).

Parameters:
  • xyzs (Iterable) – Conformers with the same atom orders.

  • rmsd_threshold (float) – The minimum RMSD to consider two conformers as distinct (i.e., if rmsd > rmsd_threshold, then two conformers are considered distinctive).

Returns:

Conformers with distinctive geometries.

Return type:

Tuple[Dict[str, tuple]]

arc.species.converter.compare_confs(xyz1: dict, xyz2: dict, rtol: float = 0.01, atol: float = 0.1, rmsd_score: bool = False) Union[float, bool][source]

Compare two Cartesian coordinates representing conformers using distance matrices.

The relative difference (rtol * abs(value in xyz2)) and the absolute difference atol are added together to compare against the absolute difference between (value in xyz1) and (value in xyz2).

Parameters:
  • xyz1 (dict) – Conformer 1.

  • xyz2 (dict) – Conformer 2.

  • rtol (float) – The relative tolerance parameter (see Notes).

  • atol (float) – The absolute tolerance parameter (see Notes).

  • rmsd_score (bool) – Whether to output a root-mean-square deviation score of the two distance matrices.

Returns:

  • If rmsd_score is False (default): Whether the two conformers have almost equal atom distances. True if they do.

  • If rmsd_score is True: The RMSD score of two distance matrices.

Return type:

Union[float, bool]

arc.species.converter.compare_zmats(z1, z2, r_tol=0.01, a_tol=2, d_tol=2, verbose=False, symmetric_torsions=None, index=1)[source]

Compare internal coordinates of two conformers of the same species. The comparison could principally be done using all dihedrals, which is information this module readily has, but this function uses Z matrices instead for better robustness (this way rings are considered as well).

Parameters:
  • z1 (dict) – Z matrix of conformer 1.

  • z2 (dict) – Z matrix of conformer 2.

  • r_tol (float, optional) – A tolerance for comparing distances (in Angstrom).

  • a_tol (float, optional) – A tolerance for comparing angles (in degrees).

  • d_tol (float, optional) – A tolerance for comparing dihedral angles (in degrees).

  • verbose (bool, optional) – Whether to print a reason for determining the zmats are different if they are, True to print.

  • symmetric_torsions (dict, optional) – Keys are tuples scan indices (0- or 1-indexed), values are internal rotation symmetry numbers (sigma). Conformers which only differ by an integer number of 360 degrees / sigma are considered identical.

  • index (int, optional) – Either 0 or 1 to specify the starting index in the keys of symmetric_torsions

Returns:

Whether the coordinates represent the same conformer within the given tolerance, True if they do.

Return type:

bool

Raises:

InputError – If xyz1 and xyz2 are of wrong type.

arc.species.converter.displace_xyz(xyz: dict, displacement: ndarray, amplitude: float = 0.25, use_weights: bool = True) Tuple[dict, dict][source]

Displace the coordinates using the displacement by the requested amplitude using atom mass weights.

Parameters:
  • xyz (dict) – The coordinates.

  • displacement (list) – The corresponding xyz displacement for each atom.

  • amplitude (float, optional) – The factor multiplication for the displacement.

  • use_weights (bool, optional) – Whether to scale displacements by the square root of the respective element mass.

Returns:

The two displaced xyz’s, one for each direction (+/-) of the weighted displacement.

Return type:

Tuple[dict, dict]

arc.species.converter.elementize(atom)[source]

Convert the atom type of an RMG Atom object into its general parent element atom type (e.g., ‘S4d’ into ‘S’).

Parameters:

atom (Atom) – The atom to process.

arc.species.converter.get_center_of_mass(xyz)[source]

Get the center of mass of xyz coordinates. Assumes arc.converter.standardize_xyz_string() was already called for xyz. Note that xyz from ESS output is usually already centered at the center of mass (to some precision).

Parameters:

xyz (dict) – The xyz coordinates.

Returns:

The center of mass coordinates.

Return type:

tuple

arc.species.converter.get_element_mass_from_xyz(xyz: dict) List[float][source]

Get a list of element masses corresponding to the given xyz considering isotopes.

Parameters:

xyz (dict) – The coordinates.

Returns:

The corresponding list of mass in amu.

Return type:

List[float]

arc.species.converter.get_most_common_isotope_for_element(element_symbol)[source]

Get the most common isotope for a given element symbol.

Parameters:

element_symbol (str) – The element symbol.

Returns:

The most common isotope number for the element.

Returns None for dummy atoms (‘X’).

Return type:

int

arc.species.converter.get_xyz_radius(xyz)[source]

Determine the largest distance from the coordinate system origin attributed to one of the atoms in 3D space.

Returns:

The radius in Angstrom.

Return type:

float

arc.species.converter.get_zmat_param_value(coords: Dict[str, tuple], indices: List[int], mol: Molecule, index: int = 0) float[source]

Generates a zmat similarly to modify_coords(), but instead of modifying it, only reports on the value of a requested parameter.

Parameters:
  • coords (dict) – Either cartesian (xyz) or internal (zmat) coordinates.

  • indices (list) – The indices to change. Specifying a list of length 2, 3, or 4 will result in a bond length, angle, or a dihedral angle parameter, respectively.

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

  • index (bool, optional) – Whether the specified atoms are 0- or 1-indexed.

Returns:

The parameter value in Angstrom or degrees.

Return type:

float

arc.species.converter.get_zmat_str_var_value(zmat_str, var)[source]

Returns the value of a zmat variable from a string-represented zmat.

Parameters:
  • zmat_str (str) – The string representation of the zmat.

  • var (str) – The variable to look for.

Returns:

The value corresponding to the var.

Return type:

float

arc.species.converter.hartree_to_si(e: float, kilo: bool = True) float[source]

Convert Hartree units into J/mol or into kJ/mol.

Parameters:
  • e (float) – Energy in Hartree.

  • kilo (bool, optional) – Whether to return kJ/mol units. True by default.

arc.species.converter.ics_to_scan_constraints(ics: list, software: Optional[str] = 'gaussian') str[source]

A helper function for converting internal coordinate (ic) info into a str block which can be read as scan constraints by ESS.

Parameters:
  • ics (list) – A list of internal coordinates (ic, stored as lists of atom indices).

  • software (str, optional) – The electronic structure software.

Returns:

A str block can be read as scan constraints by ESS.

Return type:

str

arc.species.converter.modify_coords(coords: Dict[str, tuple], indices: List[int], new_value: float, modification_type: str, mol: Optional[Molecule] = None, index: int = 0, fragments: Optional[List[List[int]]] = None) Dict[str, tuple][source]

Modify either a bond length, angle, or dihedral angle in the given coordinates. The coordinates input could either be cartesian (preferred) or internal (will be first converter to cartesian, then to internal back again since a specific zmat must be created). Internal coordinates will be used for the modification (using folding and unfolding).

Specifying an ‘atom’ modification type will only translate/rotate the atom represented by the first index if the corresponding zmat parameter is changed. Specifying a ‘group’ modification type will cause the entire group connected to the first atom to translate/rotate if the corresponding zmat parameter is changed. Specifying a ‘groups’ modification type (only valid for dihedral angles) will cause the groups connected to the first two atoms to translate/rotate if the corresponding zmat parameter is changed.

Parameters:
  • coords (dict) – Either cartesian (xyz) or internal (zmat) coordinates.

  • indices (list) – The indices to change. Specifying a list of length 2, 3, or 4 will result in changing a bond length, angle, or a dihedral angle, respectively.

  • new_value (float) – The new value to set (in Angstrom or degrees).

  • modification_type (str) – Either ‘atom’, ‘group’, or ‘groups’ (‘groups’ is only allowed for dihedral angles). Note that D ‘groups’ is a composite constraint, equivalent to calling D ‘group’ for each 1st neighboring atom in a torsion top.

  • mol (Molecule, optional) – The corresponding RMG molecule with the connectivity information. Mandatory if the modification type is ‘group’ or ‘groups’.

  • index (bool, optional) – Whether the specified atoms in indices and fragments are 0- or 1-indexed.

  • fragments (List[List[int]], optional) – Fragments represented by the species, i.e., as in a VdW well or a TS. Entries are atom index lists of all atoms in a fragment, each list represents a different fragment. indices are 0-indexed.

Raises:

InputError – If coords is not give, or if a group/s modification type is requested but mol is None, or if a ‘groups’ modification type was specified for ‘R’ or ‘A’.

Returns:

The respective cartesian (xyz) coordinates reflecting the desired modification.

Return type:

dict

arc.species.converter.molecules_from_xyz(xyz: Optional[Union[dict, str]], multiplicity: Optional[int] = None, charge: int = 0) Tuple[Optional[Molecule], Optional[Molecule]][source]

Creating RMG:Molecule objects from xyz with correct atom labeling. Based on the MolGraph.perceive_smiles method. If multiplicity is given, the returned species multiplicity will be set to it.

Parameters:
  • xyz (dict) – The ARC dict format xyz coordinates of the species.

  • multiplicity (int, optional) – The species spin multiplicity.

  • charge (int, optional) – The species net charge.

Returns: Tuple[Optional[Molecule], Optional[Molecule]]
  • The respective Molecule object with only single bonds.

  • The respective Molecule object with perceived bond orders.

arc.species.converter.order_atoms(ref_mol, mol)[source]

Order the atoms in mol by the atom order in ref_mol.

Parameters:
  • ref_mol (Molecule) – The reference Molecule object.

  • mol (Molecule) – The Molecule object to process.

Raises:
  • SanitizationError – If atoms could not be re-ordered.

  • TypeError – If mol has a wrong type.

arc.species.converter.order_atoms_in_mol_list(ref_mol, mol_list) bool[source]

Order the atoms in all molecules of mol_list by the atom order in ref_mol.

Parameters:
  • ref_mol (Molecule) – The reference Molecule object.

  • mol_list (list) – Entries are Molecule objects whose atoms will be reordered according to the reference.

Raises:

TypeError – If ref_mol or the entries in mol_list have a wrong type.

Returns:

Whether the reordering was successful, True if it was.

Return type:

bool

arc.species.converter.pybel_to_inchi(pybel_mol, has_h=True)[source]

Convert an Open Babel molecule object to InChI

Parameters:
  • pybel_mol (OBmol) – An Open Babel molecule.

  • has_h (bool) – Whether the molecule has hydrogen atoms. True if it does.

Returns:

The respective InChI representation of the molecule.

Return type:

str

arc.species.converter.rdkit_conf_from_mol(mol: Molecule, xyz: dict) tuple[source]

Generate an RDKit Conformer object from an RMG Molecule object.

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

  • xyz (dict) – The xyz coordinates of the conformer, atoms must be ordered as in mol.

Raises:

ConverterError – if xyz is of wrong type.

Returns:

  • Conformer: An RDKit Conformer object.

  • RDMol: An RDKit Molecule object.

Return type:

tuple

arc.species.converter.relocate_zmat_dummy_atoms_to_the_end(zmat_map: dict) dict[source]

Relocate all dummy atoms in a ZMat to the end of the corresponding Cartesian coordinates atom list. Only modifies the values of the ZMat map.

Parameters:

zmat_map (dict) – The ZMat map.

Returns:

The updated ZMat map.

Return type:

dict

arc.species.converter.remove_dummies(xyz)[source]

Remove dummy (‘X’) atoms from cartesian coordinates.

Parameters:

xyz (dict, str) – The cartesian coordinate, either in a dict or str format.

Returns:

The coordinates w/o dummy atoms.

Return type:

dict

Raises:

InputError – If xyz if of wrong type.

arc.species.converter.rmg_conformer_to_xyz(conformer)[source]

Convert xyz coordinates from an rmgpy.statmech.Conformer object into the ARC dict xyz style.

Notes

Only the xyz information (symbols and coordinates) will be taken from the Conformer object. Other properties such as electronic energy will not be converted.

We also assume that we can get the isotope number by rounding the mass

Parameters:

conformer (Conformer) – An rmgpy.statmech.Conformer object containing the desired xyz coordinates

Raises:

TypeError – If conformer is not an rmgpy.statmech.Conformer object

Returns:

The ARC xyz format

Return type:

dict

arc.species.converter.rmg_mol_from_inchi(inchi: str)[source]

Generate an RMG Molecule object from InChI.

Parameters:

inchi (str) – The InChI string.

Returns:

The respective RMG Molecule object.

Return type:

Molecule

arc.species.converter.s_bonds_mol_from_xyz(xyz: dict) Optional[Molecule][source]

Create a single bonded molecule from xyz using RMG’s connect_the_dots() method.

Parameters:

xyz (dict) – The xyz coordinates.

Returns: Optional[Molecule]

The respective molecule with only single bonds.

arc.species.converter.set_multiplicity(mol, multiplicity, charge, radical_map=None)[source]

Set the multiplicity and charge of a molecule. If a radical_map, which is an RMG Molecule object with the same atom order, is given, it’ll be used to set radicals (useful if bond orders aren’t known for a molecule).

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

  • multiplicity (int) –

  • charge (int) – The species net charge.

  • radical_map (Molecule, optional) – An RMG Molecule object with the same atom order to be used as a radical map.

Raises:
  • ConverterError – If radical_map is of wrong type.

  • SpeciesError – If number of radicals and multiplicity do not match or if connectivity cannot be inferred.

arc.species.converter.set_radicals_by_map(mol, radical_map)[source]

Set radicals in mol by radical_map.

Parameters:
  • mol (Molecule) – The RMG Molecule object to process.

  • radical_map (Molecule) – An RMG Molecule object with the same atom order to be used as a radical map.

Raises:

ConverterError – If atom order does not match.

arc.species.converter.set_rdkit_dihedrals(conf, rd_mol, torsion, deg_increment=None, deg_abs=None)[source]

A helper function for setting dihedral angles using RDKit. Either deg_increment or deg_abs must be specified.

Parameters:
  • conf – The RDKit conformer with the current xyz information.

  • rd_mol – The respective RDKit molecule.

  • torsion (list, tuple) – The 0-indexed atom indices of the four atoms defining the torsion.

  • deg_increment (float, optional) – The required dihedral increment in degrees.

  • deg_abs (float, optional) – The required dihedral in degrees.

Returns:

The xyz with the new dihedral, ordered according to the map.

Return type:

dict

Raises:

ConverterError – If the dihedral cannot be set.

arc.species.converter.sort_xyz_using_indices(xyz_dict: dict, indices: Optional[List[int]]) dict[source]

Sort the tuples in an xyz dict according to the given indices.

Parameters:
  • xyz_dict (dict) – The Cartesian coordinates.

  • indices (Optional[List[int]]) – Entries are 0-indices of the desired order.

Returns:

The ordered xyz.

Return type:

dict

arc.species.converter.species_to_sdf_file(species: ARCSpecies, path: str)[source]

Write an SDF file with coordinates and connectivity information.

Parameters:
  • species (ARCSpecies) – The ARCSpecies object instance.

  • path (str) – The full path to the .sdf file that will be saved.

arc.species.converter.split_str_zmat(zmat_str) Tuple[str, Optional[str]][source]

Split a string zmat into its coordinates and variables sections.

Parameters:

zmat_str (str) – The zmat.

Returns:

The coords section and the variables section if it exists, else None.

Return type:

Tuple[str, Optional[str]]

arc.species.converter.standardize_xyz_string(xyz_str, isotope_format=None)[source]

A helper function to correct xyz string format input (string to string). Usually empty lines are added by the user either in the beginning or the end, here we remove them along with other common issues.

Parameters:
  • xyz_str (str) – The string xyz format, or a Gaussian output format.

  • isotope_format (str, optional) – The format for specifying the isotope if it is not the most abundant one. By default, isotopes will not be specified. Currently the only supported option is ‘gaussian’.

Returns:

The string xyz format in standardized format.

Return type:

str

Raises:

ConverterError – If xyz_str is of wrong type.

arc.species.converter.str_to_xyz(xyz_str: str, project_directory: Optional[str] = None) dict[source]

Convert a string xyz format to the ARC dict xyz style. Note: The xyz_str argument could also direct to a file path to parse the data from. The xyz string format may have optional Gaussian-style isotope specification, e.g.:

C(Iso=13)    0.6616514836    0.4027481525   -0.4847382281
N           -0.6039793084    0.6637270105    0.0671637135
H           -1.4226865648   -0.4973210697   -0.2238712255
H           -0.4993010635    0.6531020442    1.0853092315
H           -2.2115796924   -0.4529256762    0.4144516252
H           -1.8113671395   -0.3268900681   -1.1468957003

which will also be parsed into the ARC xyz dictionary format, e.g.:

{'symbols': ('C', 'N', 'H', 'H', 'H', 'H'),
 'isotopes': (13, 14, 1, 1, 1, 1),
 'coords': ((0.6616514836, 0.4027481525, -0.4847382281),
            (-0.6039793084, 0.6637270105, 0.0671637135),
            (-1.4226865648, -0.4973210697, -0.2238712255),
            (-0.4993010635, 0.6531020442, 1.0853092315),
            (-2.2115796924, -0.4529256762, 0.4144516252),
            (-1.8113671395, -0.3268900681, -1.1468957003))}
Parameters:
  • xyz_str (str) – The string xyz format to be converted.

  • project_directory (str, optional) – The path to the project directory.

Raises:

ConverterError – If xyz_str is not a string or does not have four space-separated entries per non-empty line.

Returns: dict

The ARC xyz format.

arc.species.converter.str_to_zmat(zmat_str)[source]

Convert a string Z Matrix format to the ARC dict zmat style. Note: The zmat_str argument could also direct to a file path to parse the data from. A typical zmat string format may look like this:

C
H       1      R1
H       1      R1       2      A1
H       1      R1       2      A1       3      D1
H       1      R1       2      A1       3      D2
A1=109.4712
D1=120.0000
D2=240.0000
R1=1.0912

The resulting zmat for the above example is:

{'symbols': ('C', 'H', 'H', 'H', 'H'),
 'coords': ((None, None, None),
            ('R_1_0', None, None),
            ('R_2_1', 'A_2_1_0', None),
            ('R_3_2', 'A_3_2_0', 'D_3_2_0_1'), ('R_4_3', 'A_4_3_0', 'D_4_3_0_2')),
 'vars': {'R_1_0': 1.0912, 'R_2_1': 1.782, 'A_2_1_0': 35.2644, 'R_3_2': 1.782, 'A_3_2_0': 35.2644,
          'D_3_2_0_1': 120.0, 'R_4_3': 1.782, 'A_4_3_0': 35.2644, 'D_4_3_0_2': 240.0},
 'map': {0: 0, 1: 1, 2: 2, 3: 3, 4: 4}}
Parameters:

zmat_str (str) – The string zmat format to be converted.

Returns:

The ARC zmat format.

Return type:

dict

Raises:

ConverterError – If zmat_str is not a string or does not have enough values per line.

arc.species.converter.to_rdkit_mol(mol, remove_h=False, sanitize=True)[source]

Convert a molecular structure to an RDKit RDMol object. Uses RDKit to perform the conversion. Perceives aromaticity. Adopted from rmgpy/molecule/converter.py

Parameters:
  • mol (Molecule) – An RMG Molecule object for the conversion.

  • remove_h (bool, optional) – Whether to remove hydrogen atoms from the molecule, True to remove.

  • sanitize (bool, optional) – Whether to sanitize the RDKit molecule, True to sanitize.

Returns:

An RDKit molecule object corresponding to the input RMG Molecule object.

Return type:

RDMol

arc.species.converter.translate_to_center_of_mass(xyz)[source]

Translate coordinates to their center of mass.

Parameters:

xyz (dict) – The 3D coordinates.

Returns:

The translated coordinates.

Return type:

dict

arc.species.converter.translate_xyz(xyz_dict: dict, translation: Tuple[float, float, float]) dict[source]

Translate xyz.

Parameters:
  • xyz_dict (dict) – The ARC xyz format.

  • translation (Tuple[float, float, float]) – The x, y, z translation vector.

Returns:

The translated xyz.

Return type:

dict

arc.species.converter.update_molecule(mol, to_single_bonds=False)[source]

Updates the molecule, useful for isomorphism comparison.

Parameters:
  • mol (Molecule) – The RMG Molecule object to process.

  • to_single_bonds (bool, optional) – Whether to convert all bonds to single bonds. True to convert.

Returns:

The updated molecule.

Return type:

Molecule

arc.species.converter.xyz_file_format_to_xyz(xyz_file: str) dict[source]

Get the ARC xyz dictionary format from an XYZ file format representation.

Parameters:

xyz_file (str) – The content of an XYZ file

Raises:

ConverterError – If cannot identify the number of atoms entry, or if it is different that the actual number.

Returns: dict

The ARC dictionary xyz format.

arc.species.converter.xyz_from_data(coords, numbers=None, symbols=None, isotopes=None) dict[source]

Get the ARC xyz dictionary format from raw data. Either numbers or symbols must be specified. If isotopes isn’t specified, the most common isotopes will be assumed for all elements.

Parameters:
  • coords (tuple, list) – The xyz coordinates.

  • numbers (tuple, list, optional) – Element nuclear charge numbers.

  • symbols (tuple, list, optional) – Element symbols.

  • isotopes (tuple, list, optional) – Element isotope numbers.

Raises:

ConverterError – If neither numbers nor symbols are specified, if both are specified, or if the input lengths aren’t consistent.

Returns:

The ARC dictionary xyz format.

Return type:

dict

arc.species.converter.xyz_to_ase(xyz_dict: dict) Atoms[source]

Convert an xyz dict to an ASE Atoms object.

Parameters:

xyz_dict (dict) – The ARC xyz format.

Returns:

The corresponding ASE Atom object.

Return type:

Type[Atoms]

arc.species.converter.xyz_to_coords_and_element_numbers(xyz: dict) Tuple[list, list][source]

Convert xyz to a coords list and an atomic number list.

Parameters:

xyz (dict) – The coordinates.

Returns:

Coords and atomic numbers.

Return type:

Tuple[list, list]

arc.species.converter.xyz_to_coords_list(xyz_dict: dict) Optional[List[List[float]]][source]

Get the coords part of an xyz dict as a (mutable) list of lists (rather than a tuple of tuples).

Parameters:

xyz_dict (dict) – The ARC xyz format.

Returns: Optional[List[List[float]]]

The coordinates.

arc.species.converter.xyz_to_dmat(xyz_dict: dict) Optional[array][source]

Convert Cartesian coordinates to a distance matrix.

Parameters:

xyz_dict (dict) – The Cartesian coordinates.

Returns:

The distance matrix.

Return type:

Optional[np.array]

arc.species.converter.xyz_to_kinbot_list(xyz_dict: dict) List[Union[str, float]][source]

Get the KinBot xyz format of a single running list of: [symbol0, x0, y0, z0, symbol1, x1, y1, z1,…]

Parameters:

xyz_dict (dict) – The ARC xyz format.

Returns: List[Union[str, float]]

The respective KinBot xyz format.

arc.species.converter.xyz_to_np_array(xyz_dict: dict) Optional[ndarray][source]

Get the coords part of an xyz dict as a numpy array.

Parameters:

xyz_dict (dict) – The ARC xyz format.

Returns: Optional[np.ndarray]

The coordinates.

arc.species.converter.xyz_to_pybel_mol(xyz: dict)[source]

Convert xyz into an Open Babel molecule object.

Parameters:

xyz (dict) – ARC’s xyz dictionary format.

Returns: Optional[OBmol]

An Open Babel molecule.

arc.species.converter.xyz_to_rmg_conformer(xyz_dict: dict) Optional[Conformer][source]

Convert the Arc dict xyz style into an rmgpy.statmech.Conformer object containing these coordinates.

Notes

Only the xyz information will be supplied to the newly created Conformer object

Parameters:

xyz_dict (dict) – The ARC dict xyz style coordinates

Returns:

An rmgpy.statmech.Conformer object containing the desired xyz coordinates.

Return type:

Optional[Conformer]

arc.species.converter.xyz_to_str(xyz_dict: dict, isotope_format: Optional[str] = None) Optional[str][source]

Convert an ARC xyz dictionary format, e.g.:

{'symbols': ('C', 'N', 'H', 'H', 'H', 'H'),
 'isotopes': (13, 14, 1, 1, 1, 1),
 'coords': ((0.6616514836, 0.4027481525, -0.4847382281),
            (-0.6039793084, 0.6637270105, 0.0671637135),
            (-1.4226865648, -0.4973210697, -0.2238712255),
            (-0.4993010635, 0.6531020442, 1.0853092315),
            (-2.2115796924, -0.4529256762, 0.4144516252),
            (-1.8113671395, -0.3268900681, -1.1468957003))}

to a string xyz format with optional Gaussian-style isotope specification, e.g.:

C(Iso=13)    0.6616514836    0.4027481525   -0.4847382281
N           -0.6039793084    0.6637270105    0.0671637135
H           -1.4226865648   -0.4973210697   -0.2238712255
H           -0.4993010635    0.6531020442    1.0853092315
H           -2.2115796924   -0.4529256762    0.4144516252
H           -1.8113671395   -0.3268900681   -1.1468957003
Parameters:
  • xyz_dict (dict) – The ARC xyz format to be converted.

  • isotope_format (str, optional) – The format for specifying the isotope if it is not the most abundant one. By default, isotopes will not be specified. Currently the only supported option is ‘gaussian’.

Raises:

ConverterError – If input is not a dict or does not have all attributes.

Returns: Optional[str]

The string xyz format.

arc.species.converter.xyz_to_turbomol_format(xyz_dict: dict, charge: Optional[int] = None, unpaired: Optional[int] = None) Optional[str][source]

Get the respective Turbomole coordinates format.

$eht charge=0 unpaired=0

Args:

xyz_dict (dict): The ARC xyz format. charge (int, optional): The net charge. unpaired (int, optional): The number of unpaired electrons.

Returns:

str: The respective Turbomole coordinates.

arc.species.converter.xyz_to_x_y_z(xyz_dict: dict) Optional[Tuple[tuple, tuple, tuple]][source]

Get the X, Y, and Z coordinates separately from the ARC xyz dictionary format.

Parameters:

xyz_dict (dict) – The ARC xyz format.

Returns: Optional[Tuple[tuple, tuple, tuple]]

The X coordinates, the Y coordinates, the Z coordinates.

arc.species.converter.xyz_to_xyz_file_format(xyz_dict: dict, comment: str = '') Optional[str][source]

Get the XYZ file format representation from the ARC xyz dictionary format. This function does not consider isotopes.

Parameters:
  • xyz_dict (dict) – The ARC xyz format.

  • comment (str, optional) – A comment to be shown in the output’s 2nd line.

Raises:

ConverterError – If xyz_dict is of wrong format or comment is a multiline string.

Returns: Optional[str]

The XYZ file format.

arc.species.converter.zmat_from_xyz(xyz: Union[dict, str], mol: Optional[Molecule] = None, constraints: Optional[dict] = None, consolidate: bool = True, consolidation_tols: Optional[dict] = None, is_ts: bool = False) dict[source]

Generate a Z matrix from xyz.

Parameters:
  • xyz (Union[dict, str]) – The cartesian coordinate, either in a dict or str format.

  • mol (Molecule, optional) – The corresponding RMG Molecule with connectivity information.

  • constraints (dict, optional) – Accepted keys are: ‘R_atom’, ‘R_group’, ‘A_atom’, ‘A_group’, ‘D_atom’, ‘D_group’, or ‘D_groups’. ‘R’, ‘A’, and ‘D’ constrain distances, angles, and dihedrals, respectively. Values are lists of atom indices (0-indexed) tuples. The atom indices order matters. Specifying ‘_atom’ will cause only the last atom in the specified list values to translate/rotate if the corresponding zmat parameter is changed. Specifying ‘_group’ will cause the entire group connected to the last atom to translate/rotate if the corresponding zmat parameter is changed. Specifying ‘_groups’ (only valid for D) will cause the groups connected to the last two atoms to translate/rotate if the corresponding parameter is changed.

  • consolidate (bool, optional) – Whether to consolidate the zmat after generation, True to consolidate.

  • consolidation_tols (dict, optional) – Keys are ‘R’, ‘A’, ‘D’, values are floats representing absolute tolerance for consolidating almost equal internal coordinates.

  • is_ts (bool, optional) – Whether this is a representation of a TS. If it is not, a mol object will be generated if not given.

Raises:

InputError – If xyz if of a wrong type.

Returns:

The Z matrix.

Return type:

dict

arc.species.converter.zmat_to_str(zmat, zmat_format='gaussian', consolidate=True)[source]

Convert a zmat to a string format.

Parameters:
  • zmat (dict) – The Z Matrix to convert.

  • zmat_format (str, optional) – The requested format to output (varies by ESS). Allowed values are: ‘gaussian’, ‘qchem’, ‘molpro’, ‘orca’, ‘cfour’, or ‘psi4’. The default format is ‘gaussian’.

  • consolidate (bool) – Whether to return a consolidated zmat (geometry optimization will be more efficient).

Returns:

The string representation of the zmat in the requested format.

Return type:

str

Raises:

ConverterError – If zmat is of wrong type or missing keys, or if zmat_format is not recognized.

arc.species.converter.zmat_to_xyz(zmat, keep_dummy=False, xyz_isotopes=None)[source]

Generate the xyz dict coordinates from a zmat dict. Most common isotopes assumed, unless a reference xyz dict is given.

Parameters:
  • zmat (dict) – The zmat.

  • keep_dummy (bool) – Whether to keep dummy atoms (‘X’), True to keep, default is False.

  • xyz_isotopes (dict) – A reference xyz dictionary to take isotope information from. Must be ordered as the original mol/xyz used to create zmat.

Returns:

The xyz cartesian coordinates.

Return type:

dict