arc.species.zmat

A module for representing and manipulating Z matrices (internal coordinates).

An example for a (consolidated) zmat representation for methane:

{'symbols': ('C', 'H', 'H', 'H', 'H'),
 'coords': ((None, None, None),
            ('R_1_0', None, None),
            ('R_2|3|4_1|2|3', 'A_2|3|4_1|2|3_0|0|0', None),
            ('R_2|3|4_1|2|3', 'A_2|3|4_1|2|3_0|0|0', 'D_3_2_0_1'),
            ('R_2|3|4_1|2|3', 'A_2|3|4_1|2|3_0|0|0', 'D_4_3_0_2')),
 'vars': {'R_1_0': 1.09125,
          'D_3_2_0_1': 120.0,
          'D_4_3_0_2': 240.0,
          'R_2|3|4_1|2|3': 1.78200,
          'A_2|3|4_1|2|3_0|0|0': 35.26439},
 'map': {0: 0, 1: 1, 2: 2, 3: 3, 4: 4},
 }

Isotope information is not saved in the zmat, it exists in the xyz dict, and can be attained using the zmat map if needed.

Note that a 180 (or 0) degree for an angle in a Z matrix is not allowed, since a GIC (Generalized Internal Coordinates) optimization has no defined derivative at 180 degrees. Instead, a dummy atom ‘X’ is used. Dihedral angles may have any value.

arc.species.zmat.add_dummy_atom(zmat: dict, xyz: dict, coords: Union[list, tuple], connectivity: dict, r_atoms: list, a_atoms: list, n: int, atom_index: int) Tuple[dict, list, int, list, list, int][source]

Add a dummy atom ‘X’ to the zmat. Also updates the r_atoms and a_atoms lists for the original (non-dummy) atom.

Parameters:
  • zmat (dict) – The zmat.

  • xyz (dict) – The xyz dict.

  • coords (Union[list, tuple]) – Just the ‘coords’ part of the xyz dict.

  • connectivity (dict) – The atoms connectivity (keys are indices in the mol/xyz).

  • r_atoms (list) – The determined r_atoms.

  • a_atoms (list) – The determined a_atoms.

  • n (int) – The 0-index of the atom in the zmat to be added.

  • atom_index (int) – The 0-index of the atom in the molecule or cartesian coordinates to be added. (n and atom_index refer to the same atom, but it might have different indices in the zmat and the molecule/xyz)

Returns:

  • The zmat.

  • The coordinates (list of tuples).

  • The updated atom index in the zmat.

  • The R atom indices.

  • The A atom indices.

  • A specific atom index to be used as the last entry of the D atom indices.

Return type:

Tuple[dict, list, int, list, list, int]

arc.species.zmat.check_atom_a_constraints(atom_index: int, constraints: dict) Tuple[Optional[tuple], Optional[str]][source]

Check angle constraints for an atom. ‘A’ constraints are a list of tuples with length 3. The first atom in an A constraint is considered “constraint” for zmat generation: its angle will be relative to the second and third atoms.

Parameters:
  • atom_index (int) – The 0-indexed atom index to check.

  • constraints (dict) – The ‘R’, ‘A’, ‘D’ constraints dict. Values are lists of constraints.

Raises:

ZMatError – If the A constraint lengths do not equal three, or if the atom is constrained more than once.

Returns:

  • The atom indices to which the atom being checked is constrained. None if it is not constrained.

  • The constraint type (‘A_atom’, or ‘A_group’). None if it is not constrained.

Return type:

Tuple[Optional[tuple], Optional[str]]

arc.species.zmat.check_atom_d_constraints(atom_index: int, constraints: dict) Tuple[Optional[tuple], Optional[str]][source]

Check dihedral angle constraints for an atom. ‘D’ constraints are a list of tuples with length 4. The first atom in a D constraint is considered “constraint” for zmat generation: its dihedral angle will be relative to the second, third, and forth atoms.

Parameters:
  • atom_index (int) – The 0-indexed atom index to check.

  • constraints (dict) – The ‘R’, ‘A’, ‘D’ constraints dict. Values are lists of constraints.

Raises:

ZMatError – If the A constraint lengths do not equal three, or if the atom is constrained more than once.

Returns:

  • The atom indices to which the atom being checked is constrained. None if it is not constrained.

  • The constraint type (‘D_atom’, ‘D_group’).

Return type:

Tuple[Optional[tuple], Optional[str]]

arc.species.zmat.check_atom_r_constraints(atom_index: int, constraints: dict) Tuple[Optional[tuple], Optional[str]][source]

Check distance constraints for an atom. ‘R’ constraints are a list of tuples with length 2. The first atom in an R constraint is considered “constraint” for zmat generation: its distance will be relative to the second atom.

Parameters:
  • atom_index (int) – The 0-indexed atom index to check.

  • constraints (dict) – The ‘R’, ‘A’, ‘D’ constraints dict. Values are lists of constraints.

Raises:

ZMatError – If the R constraint lengths do not equal two, or if the atom is constrained more than once.

Returns:

  • The atom index to which the atom being checked is constrained. None if it is not constrained.

  • The constraint type (‘R_atom’, or ‘R_group’).

Return type:

Tuple[Optional[tuple], Optional[str]]

arc.species.zmat.consolidate_zmat(zmat: dict, mol: Optional[Molecule] = None, consolidation_tols: Optional[dict] = None) dict[source]

Consolidate (almost) identical vars in the zmat.

Parameters:
  • zmat (dict) – The zmat.

  • mol (Molecule, optional) – The RMG molecule, used for atom types if not None, otherwise atom symbols from the zmat are used.

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

Returns:

The consolidated zmat.

Return type:

dict

arc.species.zmat.determine_a_atoms(zmat: Dict[str, Union[dict, tuple]], coords: Union[list, tuple], connectivity: Dict[int, List[int]], r_atoms: Optional[List[int]], n: int, atom_index: int, a_constraint: Optional[Tuple[int, int]] = None, d_constraint: Optional[Tuple[int, int, int]] = None, a_constraint_type: Optional[str] = None, trivial_assignment: bool = False, fragments: Optional[List[List[int]]] = None) Optional[List[int]][source]

Determine the atoms for defining the angle A. This should be in the form: [n, r_atoms[1], <some other atom already in the zmat>]

Parameters:
  • zmat (dict) – The zmat.

  • coords (list, tuple) – Just the ‘coords’ part of the xyz dict.

  • connectivity (dict) – The atoms connectivity (keys are indices in the mol/xyz).

  • r_atoms (list) – The determined r_atoms.

  • n (int) – The 0-index of the atom in the zmat to be added.

  • atom_index (int) – The 0-index of the atom in the molecule or cartesian coordinates to be added. (n and atom_index refer to the same atom, but it might have different indices in the zmat and the molecule/xyz)

  • a_constraint (tuple, optional) – A-type constraints. The atom indices to which the atom being checked is constrained. None if it is not constrained.

  • d_constraint (tuple, optional) – D-type constraints. The atom indices to which the atom being checked is constrained. None if it is not constrained.

  • a_constraint_type (str, optional) – The A constraint type (‘A_atom’, or ‘A_group’).

  • trivial_assignment (bool, optional) – Whether to attempt assigning atoms without considering connectivity if the connectivity assignment fails.

  • 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.

Raises:

ZMatError – If the A atoms could not be determined. indices are 0-indexed.

Returns: Optional[List[int]]

The 0-indexed z-mat A atoms.

arc.species.zmat.determine_d_atoms(zmat: Dict[str, Union[dict, tuple]], xyz: Dict[str, tuple], coords: Union[list, tuple], connectivity: Dict[int, List[int]], a_atoms: Optional[List[int]], n: int, atom_index: int, d_constraint: Optional[Tuple[int, int, int]] = None, d_constraint_type: Optional[str] = None, specific_atom: Optional[int] = None, dummy: bool = False, fragments: Optional[List[List[int]]] = None) Optional[List[int]][source]

Determine the atoms for defining the dihedral angle D. This should be in the form: [n, a_atoms[1], a_atoms[2], <some other atom already in the zmat>]

Parameters:
  • zmat (dict) – The zmat.

  • xyz (dict) – The xyz dict.

  • coords (list, tuple) – Just the ‘coords’ part of the xyz dict.

  • connectivity (dict) – The atoms connectivity (keys are indices in the mol/xyz).

  • a_atoms (list) – The determined a_atoms.

  • n (int) – The 0-index of the atom in the zmat to be added.

  • atom_index (int) – The 0-index of the atom in the molecule or cartesian coordinates to be added. (n and atom_index refer to the same atom, but it might have different indices in the zmat and the molecule/xyz)

  • d_constraint (tuple, optional) – D-type constraints. The atom indices to which the atom being checked is constrained. None if it is not constrained.

  • d_constraint_type (str, optional) – The D constraint type (‘D_atom’, or ‘D_group’).

  • specific_atom (int, optional) – A 0-index of the zmat atom to be added to a_atoms to create d_atoms.

  • dummy (bool, optional) – Whether the atom being added (n) represents a dummy atom. True if it does.

  • 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.

Raises:

ZMatError – If the A atoms could not be determined. indices are 0-indexed.

Returns:

The 0-indexed z-mat D atoms.

Return type:

Optional[List[int]]

arc.species.zmat.determine_d_atoms_from_connectivity(zmat: dict, xyz: dict, coords: Union[list, tuple], connectivity: dict, a_atoms: list, atom_index: int, dummy: bool = False, allow_a_to_be_dummy: bool = False) list[source]

A helper function to determine d_atoms from the connectivity information.

Parameters:
  • zmat (dict) – The zmat.

  • xyz (dict) – The xyz dict.

  • coords (Union[list, tuple]) – Just the ‘coords’ part of the xyz dict.

  • connectivity (dict) – The atoms connectivity (keys are indices in the mol/xyz).

  • a_atoms (list) – The determined a_atoms.

  • atom_index (int) – The 0-index of the atom in the molecule or cartesian coordinates to be added. (n and atom_index refer to the same atom, but it might have different indices in the zmat and the molecule/xyz)

  • dummy (bool, optional) – Whether the atom being added (n) represents a dummy atom. True if it does.

  • allow_a_to_be_dummy (bool, optional) – Whether the last atom (‘A’) in d_atoms is allowed to by a dummy atom. True if it is, False by default.

Returns:

The d_atoms.

Return type:

list

arc.species.zmat.determine_d_atoms_without_connectivity(zmat: dict, coords: Union[list, tuple], a_atoms: list, n: int) list[source]

A helper function to determine d_atoms without connectivity information.

Parameters:
  • zmat (dict) – The zmat.

  • coords (Union[list, tuple]) – Just the ‘coords’ part of the xyz dict.

  • a_atoms (list) – The determined a_atoms.

  • n (int) – The 0-index of the atom in the zmat to be added.

Returns:

The d_atoms.

Return type:

list

arc.species.zmat.determine_r_atoms(zmat: Dict[str, Union[dict, tuple]], xyz: Dict[str, tuple], connectivity: Dict[int, List[int]], n: int, atom_index: int, r_constraint: Optional[Tuple[int]] = None, a_constraint: Optional[Tuple[int, int]] = None, d_constraint: Optional[Tuple[int, int, int]] = None, trivial_assignment: bool = False, fragments: Optional[List[List[int]]] = None) Optional[List[int]][source]

Determine the atoms for defining the distance R. This should be in the form: [n, <some other atom already in the zmat>]

Parameters:
  • zmat (dict) – The zmat.

  • xyz (dict) – The xyz dict.

  • connectivity (dict) – The atoms connectivity (keys are indices in the mol/xyz).

  • n (int) – The 0-index of the atom in the zmat to be added.

  • atom_index (int) – The 0-index of the atom in the molecule or cartesian coordinates to be added. (n and atom_index refer to the same atom, but it might have different indices in the zmat and the molecule/xyz)

  • r_constraint (tuple, optional) – R-type constraints. The atom index to which the atom being checked is constrained. None if it is not constrained.

  • a_constraint (tuple, optional) – A-type constraints. The atom indices to which the atom being checked is constrained. None if it is not constrained.

  • d_constraint (tuple, optional) – D-type constraints. The atom indices to which the atom being checked is constrained. None if it is not constrained.

  • trivial_assignment (bool, optional) – Whether to attempt assigning atoms without considering connectivity if the connectivity assignment fails.

  • 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:

ZMatError – If the R atoms could not be determined.

Returns: Optional[List[int]]

The 0-indexed z-mat R atoms.

arc.species.zmat.get_all_neighbors(mol: Molecule, atom_index: int) List[int][source]

Get atom indices of all neighbors of an atom in a molecule.

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

  • atom_index (int) – The index of the atom whose neighbors are requested.

Returns:

Atom indices of all neighbors of the requested atom.

Return type:

List[int]

arc.species.zmat.get_atom_connectivity_from_mol(mol: Molecule, atom1: Atom) List[int][source]

Get the connectivity of atom in mol. Returns heavy (non-H) atoms first.

Parameters:
  • mol (Molecule) – The molecule with connectivity information.

  • atom1 (Atom) – The atom to check connectivity for.

Returns:

0-indices of atoms in mol connected to atom.

Return type:

List[int]

arc.species.zmat.get_atom_indices_from_zmat_parameter(param: str) tuple[source]

Get the atom indices from a zmat parameter.

Examples

‘R_0_2’ –> ((0, 2),) ‘A_0_1_2’ –> ((0, 1, 2),) corresponding to angle 0-1-2 ‘D_0_1_2_4’ –> ((0, 1, 2, 4),) ‘R_0|0_3|4’ –> ((0, 3), (0, 4)) ‘A_0|0|0_1|1|1_2|3|4’ –> ((0, 1, 2), (0, 1, 3), (0, 1, 4)) corresponding to angles 0-1-2, 0-1-3, and 0-1-4 ‘D_0|0|0_1|1|1_2|3|4_5|6|9’ –> ((0, 1, 2, 5), (0, 1, 3, 6), (0, 1, 4, 9)) ‘RX_0_2’ –> ((0, 2),)

Parameters:

param (str) – The zmat parameter.

Returns:

Entries are tuples of indices, each describing R, A, or D parameters.

The tuple entries for R, A, and D types are of lengths 2, 3, and 4, respectively. The number of tuple entries depends on the number of consolidated parameters.

Return type:

tuple

arc.species.zmat.get_atom_order(xyz: Optional[Dict[str, tuple]] = None, mol: Optional[Molecule] = None, fragments: Optional[List[List[int]]] = None, constraints_dict: Optional[Dict[str, List[tuple]]] = None) List[int][source]

Get the order in which atoms should be added to the zmat from a 2D or a 3D representation.

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

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

  • 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.

  • constraints_dict (dict, optional) – A dictionary of atom constraints. The function will try to find an atom order in which all constrained atoms are after the atoms they are constraint to.

Returns:

The atom order, 0-indexed.

Return type:

List[int]

arc.species.zmat.get_atom_order_from_mol(mol: Molecule, fragment: Optional[List[int]] = None, constraints_dict: Optional[Dict[str, List[tuple]]] = None) List[int][source]

Get the order in which atoms should be added to the zmat from the 2D graph representation of the molecule.

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

  • fragment (List[int], optional) – Entries are 0-indexed atom indices to consider in the molecule. Only atoms within the fragment are considered.

  • constraints_dict (dict, optional) – A dictionary of atom constraints. The function will try to find an atom order in which all constrained atoms are after the atoms they are constraint to.

Returns:

The atom order, 0-indexed.

Return type:

List[int]

arc.species.zmat.get_atom_order_from_xyz(xyz: Dict[str, tuple], fragment: Optional[List[int]] = None) List[int][source]

Get the order in which atoms should be added to the zmat from the 3D geometry.

Parameters:
  • xyz (dict) – The 3D coordinates.

  • fragment (List[int], optional) – Entries are 0-indexed atom indices to consider in the molecule. Only atoms within the fragment are considered.

Returns:

The atom order, 0-indexed.

Return type:

List[int]

arc.species.zmat.get_connectivity(mol: Molecule) Dict[int, List[int]][source]

Get the connectivity information from the molecule object.

Parameters:

mol (Molecule) – The Molecule object.

Returns:

The connectivity information.

Keys are atom indices, values are tuples of respective edges, ordered with heavy atoms first. All indices are 0-indexed, corresponding to atom indices in mol (not in the zmat). None if xyz is given.

Return type:

Dict[int, List[int]]

arc.species.zmat.get_parameter_from_atom_indices(zmat: dict, indices: Union[list, tuple], xyz_indexed: bool = True) Union[str, tuple, list][source]

Get the zmat parameter from the atom indices. If indices are of length two, three, or four, an R, A, or D parameter is returned, respectively.

If a requested parameter represents an angle split by a dummy atom, combine the two dummy angles to get the original angle. In this case, a list of the two corresponding parameters will be returned.

Examples

[0, 2] –> ‘R_0_2’, or ‘R_0|1_2_5’, etc. [0, 2, 4] –> ‘A_0_2_4’, or a respective consolidated key [0, 2, 4, 9] –> ‘D_0_2_4_9’, or a respective consolidated key

Parameters:
  • zmat (dict) – The zmat.

  • indices (Union[list, tuple]) – Entries are 0-indices of atoms, list is of length 2, 3, or 4.

  • xyz_indexed (bool, optional) – Whether the atom indices relate to the xyz (and the zmat map will be used) or they already relate to the zmat. Default is True (relate to xyz).

Raises:
  • TypeError – If indices are of wrong type.

  • ZMatError – If indices has a wrong length, or not all indices are in the zmat map.

Returns: Union[str, tuple, list]

The corresponding zmat parameter.

arc.species.zmat.is_atom_in_new_fragment(atom_index: int, zmat: Dict[str, Union[dict, tuple]], fragments: Optional[List[List[int]]] = None, skip_atoms: Optional[List[int]] = None) bool[source]

Whether an atom is present in a new fragment that hasn’t been added to the zmat yet, and therefore atom assignment should not be done based on connectivity.

Parameters:
  • atom_index (int) – The 0-index of the atom in the molecule or cartesian coordinates to be added. (n and atom_index refer to the same atom, but it might have different indices in the zmat/molecule/xyz/fragments)

  • zmat (dict) – The zmat.

  • skip_atoms (list) – Atoms in the zmat map to ignore when checking fragments.

  • 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.

Returns:

Whether to consider connectivity for assigning atoms in z-mat variables.

Return type:

bool

arc.species.zmat.is_dummy(zmat: dict, zmat_index: int) bool[source]

Determine whether an atom in a zmat is a dummy atom by its zmat index.

Parameters:
  • zmat (dict) – The zmat with symbol and map information.

  • zmat_index (int) – The atom index in the zmat.

Raises:

ZMatError – If the index is invalid.

Returns:

Whether the atom represents a dummy atom ‘X’. True if it does.

Return type:

bool

arc.species.zmat.map_index_to_int(index: Union[int, str]) int[source]

Convert a zmat map value, e.g., 1 or ‘X15’, into an int, e.g., 1 or 15.

Parameters:

index (Union[int, str]) – The map index.

Returns: int

arc.species.zmat.order_fragments_by_constraints(fragments: List[List[int]], constraints_dict: Optional[Dict[str, List[tuple]]] = None) List[List[int]][source]

Get the order in which atoms should be added to the zmat from a 2D or a 3D representation.

Parameters:
  • fragments (List[List[int]]) – 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.

  • constraints_dict (dict, optional) – A dictionary of atom constraints. The function will try to find an atom order in which all constrained atoms are after the atoms they are constraint to.

Returns:

The ordered fragments list.

Return type:

List[List[int]]

arc.species.zmat.remove_1st_atom(zmat: dict) dict[source]

Remove the first atom of a zmat. Note: The first atom in ‘symbols’ with map key 0 is removed, it is not necessarily the first atom in the corresponding xyz with map value 0.

Parameters:

zmat (dict) – The zmat to process.

Returns:

The updated zmat.

Return type:

dict

arc.species.zmat.remove_1st_atom_references(zmat: dict) dict[source]

Remove references for the 1st atom (index 0) in the zmat in all coords and respective vars from the 5th atom (index 4) on. I.e., the 2nd, 3rd, and 4th atoms can include references to the 1st atom.

Parameters:

zmat (dict) – The zmat to process.

Returns:

The updated zmat.

Return type:

dict

arc.species.zmat.up_param(param: str, increment: Optional[int] = None, increment_list: Optional[List[int]] = None) str[source]

Increase the indices represented by a zmat parameter.

Parameters:
  • param (str) – The zmat parameter.

  • increment (int, optional) – The increment to increase by.

  • increment_list (list, optional) – Entries are individual indices to use when incrementing the param indices.

Raises:

ZMatError – If neither increment nor increment_list were specified, or if the increase resulted in a negative number.

Returns: str

The new parameter with increased indices.

arc.species.zmat.update_zmat_with_new_atom(zmat: dict, xyz: dict, coords: Union[list, tuple], n: int, atom_index: int, r_atoms: list, a_atoms: list, d_atoms: list, added_dummy: bool = False) dict[source]

Update the zmat with a new atom.

Parameters:
  • zmat (dict) – The zmat.

  • xyz (dict) – The xyz dict.

  • coords (Union[list, tuple]) – Just the ‘coords’ part of the xyz dict.

  • n (int) – The 0-index of the atom in the zmat to be added.

  • atom_index (int) – The 0-index of the atom in the molecule or cartesian coordinates to be added. (n and atom_index refer to the same atom, but it might have different indices in the zmat and the molecule/xyz)

  • r_atoms (list) – The R atom index descriptors.

  • a_atoms (list) – The A atom index descriptors.

  • d_atoms (list) – The D atom index descriptors.

  • added_dummy (bool, optional) – Whether a dummy atom was added as the last index of a_atoms.

Returns:

The updated zmat.

Return type:

dict

arc.species.zmat.xyz_to_zmat(xyz: Dict[str, tuple], mol: Optional[Molecule] = None, constraints: Optional[Dict[str, List[Tuple[int, ...]]]] = None, consolidate: bool = True, consolidation_tols: Optional[Dict[str, float]] = None, fragments: Optional[List[List[int]]] = None) Dict[str, tuple][source]

Generate a z-matrix from cartesian coordinates. The zmat is a dictionary with the following keys: - ‘symbols’: a tuple of strings representing atomic symbols; defined as a list and converted to a tuple - ‘coords’: a tuple of tuples representing internal coordinates; defined as a list and converted to a tuple - ‘vars’: a dictionary of all variables defined for the coordinates - ‘map’: a dictionary connecting atom indices in the zmat (keys) to atom indices in the mol/coords (values) This function assumes xyz has no dummy atoms. This function does not attempt to resolve constrain locks, and assumes only few non-circular constraints were given.

Parameters:
  • xyz (dict) – The xyz coordinates.

  • mol (Molecule, optional) – The corresponding RMG Molecule. If given, the bonding information will be used to generate a more meaningful zmat.

  • constraints (dict, optional) – Accepted keys are: ‘R_atom’, ‘R_group’, ‘A_atom’, ‘A_group’, ‘D_atom’, ‘D_group’. ‘R’, ‘A’, and ‘D’ are constrain distances, angles, and dihedrals, respectively. Values are lists of atom index tuples (0-indexed). The atom indices order matters. Specifying ‘_atom’ will cause only the first 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 first atom to translate/rotate if the corresponding zmat parameter is changed. Specifying ‘_groups’ (only valid for D) will cause the groups connected to the first two atoms to translate/rotate if the corresponding parameter is changed. Note that ‘D_groups’ should not be passed directly to this function, but to arc.species.converter.modify_coords() instead (it is a composite constraint).

  • 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.

  • 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:

ZMatError – If the zmat could not be generated.

Returns: Dict[str, tuple]

The z-matrix.

arc.species.zmat.zmat_to_coords(zmat: dict, keep_dummy: bool = False, skip_undefined: bool = False) Tuple[List[dict], List[str]][source]

Generate the cartesian coordinates from a zmat dict. Considers the zmat atomic map so the returned coordinates is ordered correctly. Most common isotopes assumed, if this is not the case, then isotopes should be reassigned to the xyz. This function assumes that all zmat variables relate to already defined atoms with a lower index in the zmat.

This function implements the SN-NeRF algorithm as described in: J. Parsons, J.B. Holmes, J.M Rojas, J. Tsai, C.E.M. Strauss, “Practical Conversion from Torsion Space to Cartesian Space for In Silico Protein Synthesis”, Journal of Computational Chemistry 2005, 26 (10), 1063-1068, https://doi.org/10.1002/jcc.20237

Tested in converterTest.py rather than zmatTest

Parameters:
  • zmat (dict) – The zmat.

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

  • skip_undefined (bool) – Whether to skip atoms with undefined variables, instead of raising an error. True to skip, default is False.

Raises:

ZMatError – If zmat is of wrong type or does not contain all keys.

Returns: Tuple[List[dict], List[str]]
  • The cartesian coordinates.

  • The atomic symbols corresponding to the coordinates.