Skip to content

units

units

Canonical unit conversions for force field parameters.

ForceField stores all parameters in a canonical unit system:

============  ========================
Parameter     Canonical unit
============  ========================
bond_k        kcal/(mol·Å²)
bond_eq       Å
angle_k       kcal/(mol·rad²)
angle_eq      degrees
torsion_k     kcal/mol
vdw_radius    Å
vdw_epsilon   kcal/mol
============  ========================

Energy convention: E = k·(x − x₀)² (no ½ factor), matching AMBER.

Loaders convert from format-native units to canonical on load. Savers convert from canonical back to format-native on save. Engines convert from canonical to engine-native at the boundary.

NewType wrappers provide zero-runtime-cost type safety: each unit quantity is a NewType alias of float, so mypy/pyright can catch mismatched conversions while CPython sees plain floats with no overhead.

mm3_bond_k_to_canonical

mm3_bond_k_to_canonical(k: float) -> KcalPerMolAngSq

Convert MM3 bond force constant (mdyn/Å) to canonical (kcal/mol/Ų).

Source code in q2mm/models/units.py
def mm3_bond_k_to_canonical(k: float) -> KcalPerMolAngSq:
    """Convert MM3 bond force constant (mdyn/Å) to canonical (kcal/mol/Ų)."""
    return KcalPerMolAngSq(k * MDYNA_TO_KCALMOLA2)

canonical_to_mm3_bond_k

canonical_to_mm3_bond_k(k: float) -> MdynPerAng

Convert canonical bond force constant (kcal/mol/Ų) to MM3 (mdyn/Å).

Source code in q2mm/models/units.py
def canonical_to_mm3_bond_k(k: float) -> MdynPerAng:
    """Convert canonical bond force constant (kcal/mol/Ų) to MM3 (mdyn/Å)."""
    return MdynPerAng(k * KCALMOLA2_TO_MDYNA)

mm3_angle_k_to_canonical

mm3_angle_k_to_canonical(k: float) -> KcalPerMolRadSq

Convert MM3 angle force constant (mdyn·Å/rad²) to canonical (kcal/mol/rad²).

Source code in q2mm/models/units.py
def mm3_angle_k_to_canonical(k: float) -> KcalPerMolRadSq:
    """Convert MM3 angle force constant (mdyn·Å/rad²) to canonical (kcal/mol/rad²)."""
    return KcalPerMolRadSq(k * MDYNA_RAD2_TO_KCALMOLRAD2)

canonical_to_mm3_angle_k

canonical_to_mm3_angle_k(k: float) -> MdynAngPerRadSq

Convert canonical angle force constant (kcal/mol/rad²) to MM3 (mdyn·Å/rad²).

Source code in q2mm/models/units.py
def canonical_to_mm3_angle_k(k: float) -> MdynAngPerRadSq:
    """Convert canonical angle force constant (kcal/mol/rad²) to MM3 (mdyn·Å/rad²)."""
    return MdynAngPerRadSq(k * KCALMOLRAD2_TO_MDYNA_RAD2)

canonical_to_openmm_bond_k

canonical_to_openmm_bond_k(k: float) -> KJPerMolAngSq

Convert canonical bond k (kcal/mol/Ų) → OpenMM custom-force k (kJ/mol/Ų).

Source code in q2mm/models/units.py
def canonical_to_openmm_bond_k(k: float) -> KJPerMolAngSq:
    """Convert canonical bond k (kcal/mol/Ų) → OpenMM custom-force k (kJ/mol/Ų)."""
    return KJPerMolAngSq(float(k) * KCAL_TO_KJ)

openmm_to_canonical_bond_k

openmm_to_canonical_bond_k(k: float) -> KcalPerMolAngSq

Convert OpenMM custom-force bond k (kJ/mol/Ų) → canonical (kcal/mol/Ų).

Source code in q2mm/models/units.py
def openmm_to_canonical_bond_k(k: float) -> KcalPerMolAngSq:
    """Convert OpenMM custom-force bond k (kJ/mol/Ų) → canonical (kcal/mol/Ų)."""
    return KcalPerMolAngSq(float(k) / KCAL_TO_KJ)

canonical_to_openmm_bond_k_nm

canonical_to_openmm_bond_k_nm(k: float) -> KJPerMolNmSq

Convert canonical bond k (kcal/mol/Ų) → kJ/(mol·nm²) without ½ convention.

For CustomBondForce expressions that use E = k·(r−r₀)² with r in nm. Differs from :func:canonical_to_openmm_harmonic_bond_k which includes the 2× factor for HarmonicBondForce's E = ½·k·(r−r₀)² convention.

Source code in q2mm/models/units.py
def canonical_to_openmm_bond_k_nm(k: float) -> KJPerMolNmSq:
    """Convert canonical bond k (kcal/mol/Ų) → kJ/(mol·nm²) without ½ convention.

    For CustomBondForce expressions that use E = k·(r−r₀)² with r in nm.
    Differs from :func:`canonical_to_openmm_harmonic_bond_k` which includes
    the 2× factor for HarmonicBondForce's E = ½·k·(r−r₀)² convention.
    """
    return KJPerMolNmSq(float(k) * KCAL_TO_KJ * 100.0)

canonical_to_openmm_angle_k

canonical_to_openmm_angle_k(k: float) -> KJPerMolRadSq

Convert canonical angle k (kcal/mol/rad²) → OpenMM custom-force k (kJ/mol/rad²).

Source code in q2mm/models/units.py
def canonical_to_openmm_angle_k(k: float) -> KJPerMolRadSq:
    """Convert canonical angle k (kcal/mol/rad²) → OpenMM custom-force k (kJ/mol/rad²)."""
    return KJPerMolRadSq(float(k) * KCAL_TO_KJ)

openmm_to_canonical_angle_k

openmm_to_canonical_angle_k(k: float) -> KcalPerMolRadSq

Convert OpenMM custom-force angle k (kJ/mol/rad²) → canonical (kcal/mol/rad²).

Source code in q2mm/models/units.py
def openmm_to_canonical_angle_k(k: float) -> KcalPerMolRadSq:
    """Convert OpenMM custom-force angle k (kJ/mol/rad²) → canonical (kcal/mol/rad²)."""
    return KcalPerMolRadSq(float(k) / KCAL_TO_KJ)

canonical_to_openmm_harmonic_bond_k

canonical_to_openmm_harmonic_bond_k(k: float) -> KJPerMolNmSq

Convert canonical bond k (kcal/mol/Ų) → HarmonicBondForce k (kJ/mol/nm²).

Accounts for E=k(x−x₀)² → E=½k(x−x₀)² convention and Å→nm.

Source code in q2mm/models/units.py
def canonical_to_openmm_harmonic_bond_k(k: float) -> KJPerMolNmSq:
    """Convert canonical bond k (kcal/mol/Ų) → HarmonicBondForce k (kJ/mol/nm²).

    Accounts for E=k(x−x₀)² → E=½k(x−x₀)² convention and Å→nm.
    """
    return KJPerMolNmSq(2.0 * float(k) * KCAL_TO_KJ * 100.0)

openmm_to_canonical_harmonic_bond_k

openmm_to_canonical_harmonic_bond_k(k: float) -> KcalPerMolAngSq

Convert HarmonicBondForce k (kJ/mol/nm²) → canonical bond k (kcal/mol/Ų).

Source code in q2mm/models/units.py
def openmm_to_canonical_harmonic_bond_k(k: float) -> KcalPerMolAngSq:
    """Convert HarmonicBondForce k (kJ/mol/nm²) → canonical bond k (kcal/mol/Ų)."""
    return KcalPerMolAngSq(float(k) / (2.0 * KCAL_TO_KJ * 100.0))

canonical_to_openmm_harmonic_angle_k

canonical_to_openmm_harmonic_angle_k(k: float) -> KJPerMolRadSq

Convert canonical angle k (kcal/mol/rad²) → HarmonicAngleForce k (kJ/mol/rad²).

Accounts for E=k(x−x₀)² → E=½k(x−x₀)² convention.

Source code in q2mm/models/units.py
def canonical_to_openmm_harmonic_angle_k(k: float) -> KJPerMolRadSq:
    """Convert canonical angle k (kcal/mol/rad²) → HarmonicAngleForce k (kJ/mol/rad²).

    Accounts for E=k(x−x₀)² → E=½k(x−x₀)² convention.
    """
    return KJPerMolRadSq(2.0 * float(k) * KCAL_TO_KJ)

openmm_to_canonical_harmonic_angle_k

openmm_to_canonical_harmonic_angle_k(k: float) -> KcalPerMolRadSq

Convert HarmonicAngleForce k (kJ/mol/rad²) → canonical angle k (kcal/mol/rad²).

Source code in q2mm/models/units.py
def openmm_to_canonical_harmonic_angle_k(k: float) -> KcalPerMolRadSq:
    """Convert HarmonicAngleForce k (kJ/mol/rad²) → canonical angle k (kcal/mol/rad²)."""
    return KcalPerMolRadSq(float(k) / (2.0 * KCAL_TO_KJ))

canonical_to_openmm_torsion_k

canonical_to_openmm_torsion_k(k: float) -> KJPerMol

Convert canonical torsion k (kcal/mol) → OpenMM torsion k (kJ/mol).

Source code in q2mm/models/units.py
def canonical_to_openmm_torsion_k(k: float) -> KJPerMol:
    """Convert canonical torsion k (kcal/mol) → OpenMM torsion k (kJ/mol)."""
    return KJPerMol(float(k) * KCAL_TO_KJ)

openmm_to_canonical_torsion_k

openmm_to_canonical_torsion_k(k: float) -> KcalPerMol

Convert OpenMM torsion k (kJ/mol) → canonical torsion k (kcal/mol).

Source code in q2mm/models/units.py
def openmm_to_canonical_torsion_k(k: float) -> KcalPerMol:
    """Convert OpenMM torsion k (kJ/mol) → canonical torsion k (kcal/mol)."""
    return KcalPerMol(float(k) / KCAL_TO_KJ)

canonical_to_openmm_epsilon

canonical_to_openmm_epsilon(eps: float) -> KJPerMol

Convert canonical vdW epsilon (kcal/mol) → OpenMM epsilon (kJ/mol).

Source code in q2mm/models/units.py
def canonical_to_openmm_epsilon(eps: float) -> KJPerMol:
    """Convert canonical vdW epsilon (kcal/mol) → OpenMM epsilon (kJ/mol)."""
    return KJPerMol(float(eps) * KCAL_TO_KJ)

openmm_to_canonical_epsilon

openmm_to_canonical_epsilon(eps: float) -> KcalPerMol

Convert OpenMM vdW epsilon (kJ/mol) → canonical epsilon (kcal/mol).

Source code in q2mm/models/units.py
def openmm_to_canonical_epsilon(eps: float) -> KcalPerMol:
    """Convert OpenMM vdW epsilon (kJ/mol) → canonical epsilon (kcal/mol)."""
    return KcalPerMol(float(eps) / KCAL_TO_KJ)

kj_to_kcal

kj_to_kcal(energy: float) -> KcalPerMol

Convert energy in kJ/mol → kcal/mol.

Source code in q2mm/models/units.py
def kj_to_kcal(energy: float) -> KcalPerMol:
    """Convert energy in kJ/mol → kcal/mol."""
    return KcalPerMol(float(energy) / KCAL_TO_KJ)

kcal_to_kj

kcal_to_kj(energy: float) -> KJPerMol

Convert energy in kcal/mol → kJ/mol.

Source code in q2mm/models/units.py
def kcal_to_kj(energy: float) -> KJPerMol:
    """Convert energy in kcal/mol → kJ/mol."""
    return KJPerMol(float(energy) * KCAL_TO_KJ)

ang_to_nm

ang_to_nm(length: float) -> Nanometer

Convert Angstroms to nanometers.

Source code in q2mm/models/units.py
def ang_to_nm(length: float) -> Nanometer:
    """Convert Angstroms to nanometers."""
    return Nanometer(float(length) * 0.1)

nm_to_ang

nm_to_ang(length: float) -> Angstrom

Convert nanometers to Angstroms.

Source code in q2mm/models/units.py
def nm_to_ang(length: float) -> Angstrom:
    """Convert nanometers to Angstroms."""
    return Angstrom(float(length) * 10.0)

rmin_half_to_sigma_nm

rmin_half_to_sigma_nm(radius: float) -> Nanometer

Convert Rmin/2 (Å) to LJ sigma (nm) for standard 12-6 NonbondedForce.

Source code in q2mm/models/units.py
def rmin_half_to_sigma_nm(radius: float) -> Nanometer:
    """Convert Rmin/2 (Å) to LJ sigma (nm) for standard 12-6 NonbondedForce."""
    return Nanometer(float(radius) * 2.0 / (2.0 ** (1.0 / 6.0)) * 0.1)

rmin_half_to_sigma

rmin_half_to_sigma(radius: float) -> Angstrom

Convert Rmin/2 (Å) to LJ sigma (Å) for 12-6 LJ potential.

Source code in q2mm/models/units.py
def rmin_half_to_sigma(radius: float) -> Angstrom:
    """Convert Rmin/2 (Å) to LJ sigma (Å) for 12-6 LJ potential."""
    return Angstrom(float(radius) * 2.0 / (2.0 ** (1.0 / 6.0)))

deg_to_rad

deg_to_rad(angle: float) -> Radians

Convert degrees to radians.

Source code in q2mm/models/units.py
def deg_to_rad(angle: float) -> Radians:
    """Convert degrees to radians."""
    return Radians(math.radians(float(angle)))

rad_to_deg

rad_to_deg(angle: float) -> Degrees

Convert radians to degrees.

Source code in q2mm/models/units.py
def rad_to_deg(angle: float) -> Degrees:
    """Convert radians to degrees."""
    return Degrees(math.degrees(float(angle)))

qm_to_canonical_bond_k

qm_to_canonical_bond_k(k: float) -> KcalPerMolAngSq

Convert QM bond force constant (Hartree/Bohr²) → canonical (kcal/mol/Ų).

Uses the Seminario two-step path: AU → mdyn/Å → kcal/mol/Ų.

Source code in q2mm/models/units.py
def qm_to_canonical_bond_k(k: float) -> KcalPerMolAngSq:
    """Convert QM bond force constant (Hartree/Bohr²) → canonical (kcal/mol/Ų).

    Uses the Seminario two-step path: AU → mdyn/Å → kcal/mol/Ų.
    """
    return KcalPerMolAngSq(float(k) * AU_BOND_K_TO_CANONICAL)

canonical_to_qm_bond_k

canonical_to_qm_bond_k(k: float) -> HartreePerBohrSq

Convert canonical bond k (kcal/mol/Ų) → QM (Hartree/Bohr²).

Source code in q2mm/models/units.py
def canonical_to_qm_bond_k(k: float) -> HartreePerBohrSq:
    """Convert canonical bond k (kcal/mol/Ų) → QM (Hartree/Bohr²)."""
    return HartreePerBohrSq(float(k) / AU_BOND_K_TO_CANONICAL)

qm_to_canonical_angle_k

qm_to_canonical_angle_k(k: float) -> KcalPerMolRadSq

Convert QM angle force constant (Hartree/rad²) → canonical (kcal/mol/rad²).

Uses the Seminario two-step path: AU → mdyn·Å/rad² → kcal/mol/rad².

Source code in q2mm/models/units.py
def qm_to_canonical_angle_k(k: float) -> KcalPerMolRadSq:
    """Convert QM angle force constant (Hartree/rad²) → canonical (kcal/mol/rad²).

    Uses the Seminario two-step path: AU → mdyn·Å/rad² → kcal/mol/rad².
    """
    return KcalPerMolRadSq(float(k) * AU_ANGLE_K_TO_CANONICAL)

canonical_to_qm_angle_k

canonical_to_qm_angle_k(k: float) -> HartreePerRadSq

Convert canonical angle k (kcal/mol/rad²) → QM (Hartree/rad²).

Source code in q2mm/models/units.py
def canonical_to_qm_angle_k(k: float) -> HartreePerRadSq:
    """Convert canonical angle k (kcal/mol/rad²) → QM (Hartree/rad²)."""
    return HartreePerRadSq(float(k) / AU_ANGLE_K_TO_CANONICAL)

qm_to_canonical_energy

qm_to_canonical_energy(e: float) -> KcalPerMol

Convert QM energy (Hartree) → canonical energy (kcal/mol).

Source code in q2mm/models/units.py
def qm_to_canonical_energy(e: float) -> KcalPerMol:
    """Convert QM energy (Hartree) → canonical energy (kcal/mol)."""
    return KcalPerMol(float(e) * HARTREE_TO_KCALMOL)

canonical_to_qm_energy

canonical_to_qm_energy(e: float) -> Hartree

Convert canonical energy (kcal/mol) → QM energy (Hartree).

Source code in q2mm/models/units.py
def canonical_to_qm_energy(e: float) -> Hartree:
    """Convert canonical energy (kcal/mol) → QM energy (Hartree)."""
    return Hartree(float(e) / HARTREE_TO_KCALMOL)

bohr_to_ang

bohr_to_ang(length: float) -> Angstrom

Convert Bohr to Angstrom.

Source code in q2mm/models/units.py
def bohr_to_ang(length: float) -> Angstrom:
    """Convert Bohr to Angstrom."""
    return Angstrom(float(length) * BOHR_TO_ANG)

ang_to_bohr

ang_to_bohr(length: float) -> Bohr

Convert Angstrom to Bohr.

Source code in q2mm/models/units.py
def ang_to_bohr(length: float) -> Bohr:
    """Convert Angstrom to Bohr."""
    return Bohr(float(length) / BOHR_TO_ANG)

hessian_kcalmola2_to_au

hessian_kcalmola2_to_au(k: float) -> HartreePerBohrSq

Convert Hessian element kcal/(mol·Å²) → Hartree/Bohr².

Source code in q2mm/models/units.py
def hessian_kcalmola2_to_au(k: float) -> HartreePerBohrSq:
    """Convert Hessian element kcal/(mol·Å²) → Hartree/Bohr²."""
    return HartreePerBohrSq(float(k) * KCALMOLA2_TO_HESSIAN_AU)

hessian_au_to_kcalmola2

hessian_au_to_kcalmola2(k: float) -> KcalPerMolAngSq

Convert Hessian element Hartree/Bohr² → kcal/(mol·Å²).

Source code in q2mm/models/units.py
def hessian_au_to_kcalmola2(k: float) -> KcalPerMolAngSq:
    """Convert Hessian element Hartree/Bohr² → kcal/(mol·Å²)."""
    return KcalPerMolAngSq(float(k) / KCALMOLA2_TO_HESSIAN_AU)

hessian_kjmolnm2_to_au

hessian_kjmolnm2_to_au(k: float) -> HartreePerBohrSq

Convert Hessian element kJ/(mol·nm²) → Hartree/Bohr².

Source code in q2mm/models/units.py
def hessian_kjmolnm2_to_au(k: float) -> HartreePerBohrSq:
    """Convert Hessian element kJ/(mol·nm²) → Hartree/Bohr²."""
    return HartreePerBohrSq(float(k) * KJMOLNM2_TO_HESSIAN_AU)

hessian_au_to_kjmolnm2

hessian_au_to_kjmolnm2(k: float) -> KJPerMolNmSq

Convert Hessian element Hartree/Bohr² → kJ/(mol·nm²).

Source code in q2mm/models/units.py
def hessian_au_to_kjmolnm2(k: float) -> KJPerMolNmSq:
    """Convert Hessian element Hartree/Bohr² → kJ/(mol·nm²)."""
    return KJPerMolNmSq(float(k) / KJMOLNM2_TO_HESSIAN_AU)

hessian_au_to_kjmola2

hessian_au_to_kjmola2(k: float) -> KJPerMolAngSq

Convert Hessian element Hartree/Bohr² → kJ/(mol·Å²).

Source code in q2mm/models/units.py
def hessian_au_to_kjmola2(k: float) -> KJPerMolAngSq:
    """Convert Hessian element Hartree/Bohr² → kJ/(mol·Å²)."""
    return KJPerMolAngSq(float(k) * HESSIAN_AU_TO_KJMOLA2)

hessian_kjmola2_to_au

hessian_kjmola2_to_au(k: float) -> HartreePerBohrSq

Convert Hessian element kJ/(mol·Å²) → Hartree/Bohr².

Source code in q2mm/models/units.py
def hessian_kjmola2_to_au(k: float) -> HartreePerBohrSq:
    """Convert Hessian element kJ/(mol·Å²) → Hartree/Bohr²."""
    return HartreePerBohrSq(float(k) / HESSIAN_AU_TO_KJMOLA2)