Skip to content

Molecule

molecule

Clean molecular structure representation for Q2MM.

Built on QCElemental for validated molecular data (symbols, geometry, charge, multiplicity, connectivity) with Q2MM-specific extensions (Hessian, detected bonds/angles, element-based matching).

DetectedBond dataclass

DetectedBond(atom_i: int, atom_j: int, elements: tuple[str, str], length: float, env_id: str = '', ff_row: int | None = None)

A bond detected from molecular geometry.

element_pair property

element_pair: tuple[str, str]

Sorted element pair for matching (e.g., ('C', 'F')).

DetectedAngle dataclass

DetectedAngle(atom_i: int, atom_j: int, atom_k: int, elements: tuple[str, str, str], value: float, env_id: str = '', ff_row: int | None = None)

An angle detected from molecular bonds.

element_triple property

element_triple: tuple[str, str, str]

Canonical element triple: (outer, center, outer) sorted by outer elements.

Q2MMMolecule dataclass

Q2MMMolecule(symbols: list[str], geometry: ndarray, atom_types: list[str] | None = None, charge: int = 0, multiplicity: int = 1, name: str = '', bond_tolerance: float = 1.3, hessian: ndarray | None = None, _bonds: list[DetectedBond] | None = None, _angles: list[DetectedAngle] | None = None)

Q2MM's internal molecular structure representation.

Wraps atomic symbols, coordinates, charge, and multiplicity with auto-detected bonds and angles. Optionally carries a Hessian matrix.

Can be created from XYZ files, QCElemental molecules, or raw data.

n_atoms property

n_atoms: int

Number of atoms in the molecule.

bonds property

bonds: list[DetectedBond]

Auto-detected bonds from covalent radii.

angles property

angles: list[DetectedAngle]

Auto-detected angles from bonds.

__post_init__

__post_init__()

Validate atom_types length and normalize geometry to float.

Source code in q2mm/models/molecule.py
def __post_init__(self):
    """Validate atom_types length and normalize geometry to float."""
    self.symbols = [str(symbol) for symbol in self.symbols]
    if self.atom_types is None:
        self.atom_types = list(self.symbols)
    else:
        self.atom_types = [str(atom_type) for atom_type in self.atom_types]
    if len(self.atom_types) != len(self.symbols):
        raise ValueError("atom_types must have the same length as symbols.")
    self.geometry = np.asarray(self.geometry, dtype=float)

from_xyz classmethod

from_xyz(path: str | Path, charge: int = 0, multiplicity: int = 1, name: str = '', bond_tolerance: float = 1.3) -> Q2MMMolecule

Load from XYZ file.

Parameters:

Name Type Description Default
bond_tolerance float

Multiplier on sum of covalent radii for bond detection. Use 1.3 for ground states, 1.4-1.5 for transition states with partially formed/broken bonds.

1.3
Source code in q2mm/models/molecule.py
@classmethod
def from_xyz(
    cls, path: str | Path, charge: int = 0, multiplicity: int = 1, name: str = "", bond_tolerance: float = 1.3
) -> Q2MMMolecule:
    """Load from XYZ file.

    Args:
        bond_tolerance: Multiplier on sum of covalent radii for bond detection.
                       Use 1.3 for ground states, 1.4-1.5 for transition states
                       with partially formed/broken bonds.
    """
    path = Path(path)
    with open(path) as f:
        lines = f.readlines()
    n = int(lines[0].strip())
    symbols = []
    coords = []
    for line in lines[2 : 2 + n]:
        parts = line.split()
        symbols.append(parts[0])
        coords.append([float(x) for x in parts[1:4]])
    return cls(
        symbols=symbols,
        atom_types=list(symbols),
        geometry=np.array(coords),
        charge=charge,
        multiplicity=multiplicity,
        name=name or path.stem,
        bond_tolerance=bond_tolerance,
    )

from_structure classmethod

from_structure(structure, *, charge: int = 0, multiplicity: int = 1, name: str = '', bond_tolerance: float = 1.3, hessian: ndarray | None = None) -> Q2MMMolecule

Create from a legacy Structure while preserving bond/angle metadata.

Source code in q2mm/models/molecule.py
@classmethod
def from_structure(
    cls,
    structure,
    *,
    charge: int = 0,
    multiplicity: int = 1,
    name: str = "",
    bond_tolerance: float = 1.3,
    hessian: np.ndarray | None = None,
) -> Q2MMMolecule:
    """Create from a legacy Structure while preserving bond/angle metadata."""
    symbols = []
    atom_types = []
    coords = []
    for atom in structure.atoms:
        atom_label = atom.atom_type_name or atom.element or ""
        symbols.append(_extract_element(atom_label))
        atom_types.append(atom_label.strip() or _extract_element(atom_label))
        coords.append(atom.coords)

    bonds = []
    for bond in structure.bonds:
        atoms = structure.get_atoms_in_DOF(bond)
        dof_atom_types = [atom.atom_type_name or atom.element or "" for atom in atoms]
        elements = tuple(_extract_element(atom_type) for atom_type in dof_atom_types[:2])
        length = bond.value
        if length is None:
            length = np.linalg.norm(atoms[0].coords - atoms[1].coords)
        bonds.append(
            DetectedBond(
                atom_i=bond.atom_nums[0] - 1,
                atom_j=bond.atom_nums[1] - 1,
                elements=elements,
                length=float(length),
                env_id=canonicalize_bond_env_id(dof_atom_types),
                ff_row=bond.ff_row,
            )
        )

    angles = []
    for angle in structure.angles:
        atoms = structure.get_atoms_in_DOF(angle)
        dof_atom_types = [atom.atom_type_name or atom.element or "" for atom in atoms]
        elements = tuple(_extract_element(atom_type) for atom_type in dof_atom_types[:3])
        angle_value = angle.value
        if angle_value is None:
            v1 = atoms[0].coords - atoms[1].coords
            v2 = atoms[2].coords - atoms[1].coords
            cos_a = np.dot(v1, v2) / (np.linalg.norm(v1) * np.linalg.norm(v2))
            angle_value = np.degrees(np.arccos(np.clip(cos_a, -1, 1)))
        angles.append(
            DetectedAngle(
                atom_i=angle.atom_nums[0] - 1,
                atom_j=angle.atom_nums[1] - 1,
                atom_k=angle.atom_nums[2] - 1,
                elements=elements,
                value=float(angle_value),
                env_id=canonicalize_angle_env_id(dof_atom_types),
                ff_row=angle.ff_row,
            )
        )

    return cls(
        symbols=symbols,
        atom_types=atom_types,
        geometry=np.array(coords, dtype=float),
        charge=charge,
        multiplicity=multiplicity,
        name=name or structure.origin_name,
        bond_tolerance=bond_tolerance,
        hessian=hessian,
        _bonds=bonds,
        _angles=angles,
    )

from_qcel classmethod

from_qcel(mol: Molecule, name: str = '') -> Q2MMMolecule

Create from a QCElemental Molecule object.

Source code in q2mm/models/molecule.py
@classmethod
def from_qcel(cls, mol: qcel.models.Molecule, name: str = "") -> Q2MMMolecule:
    """Create from a QCElemental Molecule object."""
    if not _HAS_QCEL:
        raise ImportError("qcelemental required: pip install qcelemental")
    coords_bohr = np.array(mol.geometry).reshape(-1, 3)
    coords_ang = coords_bohr * qcel.constants.bohr2angstroms
    return cls(
        symbols=list(mol.symbols),
        atom_types=list(mol.symbols),
        geometry=coords_ang,
        charge=mol.molecular_charge,
        multiplicity=mol.molecular_multiplicity,
        name=name,
    )

to_qcel

to_qcel() -> Molecule

Convert to QCElemental Molecule.

Source code in q2mm/models/molecule.py
def to_qcel(self) -> qcel.models.Molecule:
    """Convert to QCElemental Molecule."""
    if not _HAS_QCEL:
        raise ImportError("qcelemental required: pip install qcelemental")
    coords_bohr = self.geometry / qcel.constants.bohr2angstroms
    conn = [(b.atom_i, b.atom_j, 1) for b in self.bonds]
    kwargs = {
        "symbols": self.symbols,
        "geometry": coords_bohr.flatten().tolist(),
        "molecular_charge": self.charge,
        "molecular_multiplicity": self.multiplicity,
    }
    if conn:
        kwargs["connectivity"] = conn
    return qcel.models.Molecule(**kwargs)

with_hessian

with_hessian(hessian: ndarray) -> Q2MMMolecule

Return a copy with Hessian attached.

Source code in q2mm/models/molecule.py
def with_hessian(self, hessian: np.ndarray) -> Q2MMMolecule:
    """Return a copy with Hessian attached."""
    return Q2MMMolecule(
        symbols=self.symbols,
        atom_types=list(self.atom_types),
        geometry=self.geometry.copy(),
        charge=self.charge,
        multiplicity=self.multiplicity,
        name=self.name,
        bond_tolerance=self.bond_tolerance,
        hessian=hessian,
        _bonds=copy.deepcopy(self._bonds) if self._bonds is not None else None,
        _angles=copy.deepcopy(self._angles) if self._angles is not None else None,
    )