torchff#
Bond#
Bond energy and force terms (harmonic and Amoeba-style) for force field calculations.
- torchff.bond.compute_harmonic_bond_energy(coords: Tensor, bonds: Tensor, b0: Tensor, k: Tensor) Tensor[source]#
Compute harmonic bond energies via custom CUDA/CPU ops.
Energy per bond: \(E = (1/2) k (r - r_0)^2\) with \(r\) the bond length.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
bonds (torch.Tensor) – Shape (M, 2), integer indices; each row is (i, j).
b0 (torch.Tensor) – Shape (M,), equilibrium bond lengths.
k (torch.Tensor) – Shape (M,), force constants.
- Returns:
Scalar total harmonic bond energy.
- Return type:
- torchff.bond.compute_harmonic_bond_energy_ref(coords, bonds, b0, k)[source]#
Reference implementation of harmonic bond energy (PyTorch only).
Same formula as
compute_harmonic_bond_energy(); used when custom ops are disabled.
- class torchff.bond.HarmonicBond(use_customized_ops: bool = False)[source]#
Bases:
ModuleHarmonic bond energy module: \(E = (1/2) k (r - r_0)^2\).
Dispatches to custom ops or a PyTorch reference implementation based on
use_customized_ops.- __init__(use_customized_ops: bool = False)[source]#
- Parameters:
use_customized_ops (bool, optional) – If True, use custom CUDA/CPU kernels; otherwise use reference implementation.
- forward(coords, bonds, b0, k)[source]#
Compute total harmonic bond energy.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
bonds (torch.Tensor) – Shape (M, 2), bond indices.
b0 (torch.Tensor) – Shape (M,), equilibrium lengths.
k (torch.Tensor) – Shape (M,), force constants.
- Returns:
Scalar total energy.
- Return type:
- torchff.bond.compute_harmonic_bond_forces(coords: Tensor, bonds: Tensor, b0: Tensor, k: Tensor, forces: Tensor)[source]#
Compute harmonic bond forces in-place (no autograd).
Adds bond force contributions into
forces. Intended for fast MD; backward is not supported.- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
bonds (torch.Tensor) – Shape (M, 2), bond indices.
b0 (torch.Tensor) – Shape (M,), equilibrium lengths.
k (torch.Tensor) – Shape (M,), force constants.
forces (torch.Tensor) – Shape (N, 3), modified in-place with bond forces.
- torchff.bond.compute_amoeba_bond_energy(coords: Tensor, bonds: Tensor, b0: Tensor, k: Tensor, cubic: float = -2.55, quartic: float = 3.793125) Tensor[source]#
Compute Amoeba-style bond energies via custom ops.
\(E = k (b - b_0)^2 [1 + c_3 (b - b_0) + c_4 (b - b_0)^2]\) with \(b\) the bond length and \(c_3, c_4\) the cubic and quartic coefficients.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
bonds (torch.Tensor) – Shape (M, 2), bond indices.
b0 (torch.Tensor) – Shape (M,), equilibrium lengths.
k (torch.Tensor) – Shape (M,), force constants.
cubic (float, optional) – Cubic coefficient in the polynomial (default -2.55).
quartic (float, optional) – Quartic coefficient (default 3.793125).
- Returns:
Scalar total Amoeba bond energy.
- Return type:
- torchff.bond.compute_amoeba_bond_energy_ref(coords, bonds, b0, k, cubic: float = -2.55, quartic: float = 3.793125)[source]#
Reference implementation of Amoeba bond energy (PyTorch only).
Same formula as
compute_amoeba_bond_energy(); used when custom ops are disabled.
- class torchff.bond.AmoebaBond(use_customized_ops: bool = False)[source]#
Bases:
ModuleAmoeba-style bond energy with cubic and quartic correction terms.
Dispatches to custom ops or reference based on
use_customized_ops. SeeHarmonicBondfor constructor and dispatch behavior.- __init__(use_customized_ops: bool = False)[source]#
Same as
HarmonicBond.__init__.
- forward(coords, bonds, b0, k, cubic: float = -2.55, quartic: float = 3.793125)[source]#
Compute total Amoeba bond energy.
- Parameters:
coords – Same as
HarmonicBond.forward().bonds – Same as
HarmonicBond.forward().b0 – Same as
HarmonicBond.forward().k – Same as
HarmonicBond.forward().cubic (float, optional) – Cubic coefficient (default -2.55).
quartic (float, optional) – Quartic coefficient (default 3.793125).
- Returns:
Scalar total energy.
- Return type:
- torchff.bond.compute_morse_bond_energy(coords: Tensor, bonds: Tensor, b0: Tensor, k: Tensor, d: Tensor) Tensor[source]#
Compute Morse bond energies via custom CUDA/CPU ops.
Energy per bond: \(E = D [1 - e^{-\beta (r - r_0)}]^2\) with \(\beta = \sqrt{k/(2D)}\), \(r\) the bond length, \(r_0\) the equilibrium length, \(k\) the stiffness, and \(D\) the well depth.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
bonds (torch.Tensor) – Shape (M, 2), integer indices; each row is (i, j).
b0 (torch.Tensor) – Shape (M,), equilibrium bond lengths \(r_0\).
k (torch.Tensor) – Shape (M,), stiffness parameters \(k\).
d (torch.Tensor) – Shape (M,), well depth \(D\) (must be positive).
- Returns:
Scalar total Morse bond energy.
- Return type:
- torchff.bond.compute_morse_bond_energy_ref(coords, bonds, b0, k, d)[source]#
Reference implementation of Morse bond energy (PyTorch only).
Same formula as
compute_morse_bond_energy(); used when custom ops are disabled.
- class torchff.bond.MorseBond(use_customized_ops: bool = False)[source]#
Bases:
ModuleMorse bond energy module.
Dispatches to custom ops or a PyTorch reference implementation based on
use_customized_ops.- forward(coords, bonds, b0, k, d)[source]#
Compute total Morse bond energy.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
bonds (torch.Tensor) – Shape (M, 2), bond indices.
b0 (torch.Tensor) – Shape (M,), equilibrium lengths.
k (torch.Tensor) – Shape (M,), stiffness parameters.
d (torch.Tensor) – Shape (M,), well depths.
- Returns:
Scalar total energy.
- Return type:
- torchff.bond.compute_morse_bond_forces(coords: Tensor, bonds: Tensor, b0: Tensor, k: Tensor, d: Tensor, forces: Tensor)[source]#
Compute Morse bond forces in-place (no autograd).
Adds bond force contributions into
forces.- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
bonds (torch.Tensor) – Shape (M, 2), bond indices.
b0 (torch.Tensor) – Shape (M,), equilibrium lengths.
k (torch.Tensor) – Shape (M,), stiffness parameters.
d (torch.Tensor) – Shape (M,), well depths.
forces (torch.Tensor) – Shape (N, 3), modified in-place with bond forces.
Angle#
Angle energy and force terms (harmonic and Amoeba-style) for force field calculations.
- torchff.angle.compute_angles_ref(coords: Tensor, angles: Tensor) Tensor[source]#
Reference implementation: angle (theta) at the central atom for each triple (PyTorch only).
For each angle (i, j, k), computes the angle at j between vectors (coords[i] - coords[j]) and (coords[k] - coords[j]).
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
angles (torch.Tensor) – Shape (M, 3), integer indices; each row is (i, j, k).
- Returns:
Shape (M,), angle in radians for each triple.
- Return type:
- torchff.angle.compute_angles(coords: Tensor, angles: Tensor) Tensor[source]#
Compute angle values via custom CUDA/CPU ops.
Same geometry as
compute_angles_ref().- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
angles (torch.Tensor) – Shape (M, 3), integer indices; each row is (i, j, k).
- Returns:
Shape (M,), angle in radians for each triple.
- Return type:
- class torchff.angle.Angles(use_customized_ops: bool = False)[source]#
Bases:
ModulePlane angle values at central atoms (radians), matching
compute_angles_ref()geometry.Dispatches to custom ops or reference based on
use_customized_ops.
- torchff.angle.compute_harmonic_angle_energy(coords: Tensor, angles: Tensor, theta0: Tensor, k: Tensor) Tensor[source]#
Compute harmonic angle energies via custom CUDA/CPU ops.
Energy per angle: \(E = (1/2) k (\theta - \theta_0)^2\) with \(\theta\) in radians.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
angles (torch.Tensor) – Shape (M, 3), integer indices; each row is (i, j, k), j the central atom.
theta0 (torch.Tensor) – Shape (M,), equilibrium angles in radians.
k (torch.Tensor) – Shape (M,), force constants.
- Returns:
Scalar total harmonic angle energy.
- Return type:
- torchff.angle.compute_harmonic_angle_energy_ref(coords, angles, theta0, k)[source]#
Reference implementation of harmonic angle energy (PyTorch only).
Uses
compute_angles_ref()then same formula ascompute_harmonic_angle_energy(); used when custom ops are disabled.
- class torchff.angle.HarmonicAngle(use_customized_ops: bool = False)[source]#
Bases:
ModuleHarmonic angle energy module: \(E = (1/2) k (\theta - \theta_0)^2\).
Dispatches to custom ops or a PyTorch reference implementation based on
use_customized_ops.- __init__(use_customized_ops: bool = False)[source]#
- Parameters:
use_customized_ops (bool, optional) – If True, use custom CUDA/CPU kernels; otherwise use reference implementation.
- forward(coords, angles, theta0, k)[source]#
Compute total harmonic angle energy.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
angles (torch.Tensor) – Shape (M, 3), angle indices (i, j, k).
theta0 (torch.Tensor) – Shape (M,), equilibrium angles in radians.
k (torch.Tensor) – Shape (M,), force constants.
- Returns:
Scalar total energy.
- Return type:
- torchff.angle.compute_harmonic_angle_forces(coords: Tensor, angles: Tensor, theta0: Tensor, k: Tensor, forces: Tensor) Tensor[source]#
Compute harmonic angle forces in-place (no autograd).
Adds angle force contributions into
forces. Intended for fast MD; backward is not supported.- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
angles (torch.Tensor) – Shape (M, 3), angle indices.
theta0 (torch.Tensor) – Shape (M,), equilibrium angles in radians.
k (torch.Tensor) – Shape (M,), force constants.
forces (torch.Tensor) – Shape (N, 3), modified in-place with angle forces.
- Returns:
Reference to
forces.- Return type:
- torchff.angle.compute_amoeba_angle_energy(coords: Tensor, angles: Tensor, theta0: Tensor, k: Tensor, cubic: float = -0.014, quartic: float = 5.6e-05, pentic: float = -7e-07, sextic: float = 2.2e-08) Tensor[source]#
Compute Amoeba-style angle energies via custom ops.
\(E = k (\theta - \theta_0)^2 [1 + c_3 \Delta\theta + c_4 \Delta\theta^2 + c_5 \Delta\theta^3 + c_6 \Delta\theta^4]\) with \(\Delta\theta = \theta - \theta_0\) in radians.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
angles (torch.Tensor) – Shape (M, 3), angle indices (i, j, k).
theta0 (torch.Tensor) – Shape (M,), equilibrium angles in radians.
k (torch.Tensor) – Shape (M,), force constants.
cubic (float, optional) – Polynomial coefficients (defaults: -0.014, 5.6e-5, -7.0e-7, 2.2e-8).
quartic (float, optional) – Polynomial coefficients (defaults: -0.014, 5.6e-5, -7.0e-7, 2.2e-8).
pentic (float, optional) – Polynomial coefficients (defaults: -0.014, 5.6e-5, -7.0e-7, 2.2e-8).
sextic (float, optional) – Polynomial coefficients (defaults: -0.014, 5.6e-5, -7.0e-7, 2.2e-8).
- Returns:
Scalar total Amoeba angle energy.
- Return type:
- torchff.angle.compute_amoeba_angle_energy_ref(coords, angles, theta0, k, cubic: float = -0.014, quartic: float = 5.6e-05, pentic: float = -7e-07, sextic: float = 2.2e-08)[source]#
Reference implementation of Amoeba angle energy (PyTorch only).
Uses
compute_angles_ref(); same formula ascompute_amoeba_angle_energy(); used when custom ops are disabled.
- class torchff.angle.AmoebaAngle(use_customized_ops: bool = False)[source]#
Bases:
ModuleAmoeba-style angle energy with cubic through sextic correction terms.
Dispatches to custom ops or reference based on
use_customized_ops. SeeHarmonicAnglefor constructor and dispatch behavior.- __init__(use_customized_ops: bool = False)[source]#
Same as
HarmonicAngle.__init__.
- forward(coords, angles, theta0, k, cubic: float = -0.014, quartic: float = 5.6e-05, pentic: float = -7e-07, sextic: float = 2.2e-08)[source]#
Compute total Amoeba angle energy.
- Parameters:
coords – Same as
HarmonicAngle.forward().angles – Same as
HarmonicAngle.forward().theta0 – Same as
HarmonicAngle.forward().k – Same as
HarmonicAngle.forward().cubic (float, optional) – Polynomial coefficients (defaults: -0.014, 5.6e-5, -7.0e-7, 2.2e-8).
quartic (float, optional) – Polynomial coefficients (defaults: -0.014, 5.6e-5, -7.0e-7, 2.2e-8).
pentic (float, optional) – Polynomial coefficients (defaults: -0.014, 5.6e-5, -7.0e-7, 2.2e-8).
sextic (float, optional) – Polynomial coefficients (defaults: -0.014, 5.6e-5, -7.0e-7, 2.2e-8).
- Returns:
Scalar total energy.
- Return type:
- torchff.angle.compute_cosine_angle_energy(coords: Tensor, angles: Tensor, theta0: Tensor, k: Tensor) Tensor[source]#
Cosine-style angle energy via custom ops.
Energy per angle: \(E = (1/2) k (\cos\theta - \cos\theta_0)^2\).
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
angles (torch.Tensor) – Shape (M, 3), angle indices (i, j, k).
theta0 (torch.Tensor) – Shape (M,), equilibrium angles in radians.
k (torch.Tensor) – Shape (M,), force constants.
- Returns:
Scalar total cosine angle energy.
- Return type:
- torchff.angle.compute_cosine_angle_energy_ref(coords, angles, theta0, k)[source]#
Reference implementation of cosine angle energy (PyTorch only).
Same formula as
compute_cosine_angle_energy().
- class torchff.angle.CosineAngle(use_customized_ops: bool = False)[source]#
Bases:
ModuleCosine angle energy: \(E = (1/2) k (\cos\theta - \cos\theta_0)^2\).
Dispatches to custom ops or reference based on
use_customized_ops.
- torchff.angle.compute_cosine_angle_forces(coords: Tensor, angles: Tensor, theta0: Tensor, k: Tensor, forces: Tensor) Tensor[source]#
Cosine angle forces in-place (no autograd).
- Parameters:
coords – Same as
compute_cosine_angle_energy().angles – Same as
compute_cosine_angle_energy().theta0 – Same as
compute_cosine_angle_energy().k – Same as
compute_cosine_angle_energy().forces (torch.Tensor) – Shape (N, 3), modified in-place.
Torsion#
Torsion (dihedral) angles, harmonic torsion energy, and periodic torsion energies.
- torchff.torsion.compute_torsion_ref(coords: Tensor, torsions: Tensor) Tensor[source]#
Reference implementation of torsion (dihedral) angles using PyTorch ops.
For each torsion (i, j, k, l), computes the dihedral angle between the planes formed by (i, j, k) and (j, k, l). The convention matches the CUDA implementation used by the custom operator.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
torsions (torch.Tensor) – Shape (M, 4), integer indices; each row is (i, j, k, l).
- Returns:
Shape (M,), torsion angles in radians.
- Return type:
- torchff.torsion.compute_torsion(coords: Tensor, torsions: Tensor) Tensor[source]#
Compute torsion (dihedral) angles via custom CUDA/CPU ops.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
torsions (torch.Tensor) – Shape (M, 4), integer indices; each row is (i, j, k, l).
- Returns:
Shape (M,), torsion angles in radians.
- Return type:
- class torchff.torsion.Torsion(use_customized_ops: bool = False)[source]#
Bases:
ModuleTorsion (dihedral) angle module returning raw angle values in radians.
Dispatches to custom ops or a PyTorch reference implementation based on
use_customized_ops.- __init__(use_customized_ops: bool = False)[source]#
- Parameters:
use_customized_ops (bool, optional) – If True, use custom CUDA/CPU kernels; otherwise use reference implementation.
- forward(coords: Tensor, torsions: Tensor) Tensor[source]#
Compute torsion (dihedral) angles.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
torsions (torch.Tensor) – Shape (M, 4), torsion indices (i, j, k, l).
- Returns:
Shape (M,), torsion angles in radians.
- Return type:
- torchff.torsion.compute_periodic_torsion_energy(coords: Tensor, torsions: Tensor, fc: Tensor, periodicity: Tensor, phase: Tensor) Tensor[source]#
Compute periodic torsion energies
- torchff.torsion.compute_periodic_torsion_energy_ref(coords, torsions, fc, periodicity, phase)[source]#
- torchff.torsion.compute_harmonic_torsion_energy_ref(coords: Tensor, torsions: Tensor, k: Tensor, torsion_ref: Tensor | None = None) Tensor[source]#
Reference harmonic torsion energy using PyTorch dihedral angles.
Total energy is
sum_i (k_i / 2) * (phi_i - phi0_i)^2whentorsion_ref(phi0) is given; iftorsion_refis None, usessum_i (k_i / 2) * phi_i^2.- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
torsions (torch.Tensor) – Shape (M, 4), integer indices (i, j, k, l) per torsion.
k (torch.Tensor) – Shape (M,), force constants.
torsion_ref (torch.Tensor, optional) – Shape (M,), equilibrium torsion angles in radians. If None, the equilibrium angle is taken as zero for every term.
- Returns:
Scalar total energy.
- Return type:
- class torchff.torsion.HarmonicTorsion(use_customized_ops: bool = False)[source]#
Bases:
ModuleHarmonic torsion energy: sum_i (k_i/2) (phi_i - phi0_i)^2.
Dihedral angles phi_i come from
Torsion(custom or reference ops). The energy reduction is pure PyTorch; there is no separate fused CUDA kernel.- __init__(use_customized_ops: bool = False)[source]#
- Parameters:
use_customized_ops (bool, optional) – If True, dihedral angles use
compute_torsion(); otherwisecompute_torsion_ref().
- forward(coords: Tensor, torsions: Tensor, k: Tensor, torsion_ref: Tensor | None = None) Tensor[source]#
Compute total harmonic torsion energy.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
torsions (torch.Tensor) – Shape (M, 4), torsion indices (i, j, k, l).
k (torch.Tensor) – Shape (M,), force constants.
torsion_ref (torch.Tensor, optional) – Shape (M,), equilibrium angles in radians. If None, uses
sum (k/2) phi^2.
- Returns:
Scalar total energy.
- Return type:
Coulomb#
- torchff.coulomb.compute_coulomb_energy(coords: Tensor, pairs: Tensor, box: Tensor, charges: Tensor, coulomb_constant: float, cutoff: float, do_shift: bool = True) Tensor[source]#
Compute Coulomb energies using customized CUDA/C++ ops
Van der Waals#
Van der Waals (vdW) pair energies: Lennard-Jones 12-6 and AMOEBA buffered 14-7.
- torchff.vdw.compute_vdw_14_7_energy(coords: Tensor, pairs: Tensor, box: Tensor, sigma: Tensor, epsilon: Tensor, cutoff: float, atom_types: Tensor | None = None) Tensor[source]#
Compute AMOEBA buffered 14-7 vdW pair energies via custom CUDA/C++ ops.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
pairs (torch.Tensor) – Shape (P, 2), integer indices (i, j) of interacting pairs.
box (torch.Tensor) – Shape (3, 3) or broadcastable, periodic box (same convention as
torchff.pbc).sigma (torch.Tensor) – Per-pair or type-pair \(\sigma\) (see
atom_types).epsilon (torch.Tensor) – Per-pair or type-pair \(\epsilon\) (see
atom_types).cutoff (float) – Distance cutoff; interactions beyond cutoff are excluded by the kernel.
atom_types (torch.Tensor, optional) – If provided, used by the backend for type-based indexing together with
sigmaandepsilon.
- Returns:
Scalar total vdW energy for the buffered 14-7 potential.
- Return type:
- torchff.vdw.compute_lennard_jones_energy(coords: Tensor, pairs: Tensor, box: Tensor, sigma: Tensor, epsilon: Tensor, cutoff: float, atom_types: Tensor | None = None) Tensor[source]#
Compute Lennard-Jones 12-6 vdW pair energies via custom CUDA/C++ ops.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
pairs (torch.Tensor) – Shape (P, 2), integer indices (i, j) of interacting pairs.
box (torch.Tensor) – Shape (3, 3) or broadcastable, periodic box (same convention as
torchff.pbc).sigma (torch.Tensor) – Per-pair or type-pair \(\sigma\) (see
atom_types).epsilon (torch.Tensor) – Per-pair or type-pair \(\epsilon\) (see
atom_types).cutoff (float) – Distance cutoff; interactions beyond cutoff are excluded by the kernel.
atom_types (torch.Tensor, optional) – If provided, used by the backend for type-based indexing together with
sigmaandepsilon.
- Returns:
Scalar total Lennard-Jones energy.
- Return type:
- torchff.vdw.compute_lennard_jones_energy_ref(r_ij, sigma_ij, epsilon_ij, sum=True)[source]#
Reference Lennard-Jones 12-6 pair energy in PyTorch.
Per pair: \(E_{ij} = 4 \epsilon_{ij} \left[ (\sigma_{ij}/r_{ij})^{12} - (\sigma_{ij}/r_{ij})^6 \right]\).
- Parameters:
r_ij (torch.Tensor) – Pair distances, shape (P,) or broadcastable.
sigma_ij (torch.Tensor) – \(\sigma\) for each pair, same shape as
r_ij(after broadcast).epsilon_ij (torch.Tensor) – \(\epsilon\) for each pair, same shape as
r_ij(after broadcast).sum (bool, optional) – If True (default), return the sum over pairs; otherwise return per-pair energies.
- Returns:
Scalar total energy if
sumis True, else shape (P,) per-pair energies.- Return type:
- torchff.vdw.compute_vdw_14_7_energy_ref(r_ij, sigma_ij, epsilon_ij, sum=True)[source]#
Reference AMOEBA buffered 14-7 vdW pair energy in PyTorch.
With \(\rho = r_{ij} / \sigma_{ij}\), \(E_{ij} = \epsilon_{ij} \left( \frac{1.07}{\rho + 0.07} \right)^7 \left( \frac{1.12}{\rho^7 + 0.12} - 2 \right)\).
- Parameters:
r_ij (torch.Tensor) – Pair distances, shape (P,) or broadcastable.
sigma_ij (torch.Tensor) – \(\sigma\) for each pair, same shape as
r_ij(after broadcast).epsilon_ij (torch.Tensor) – \(\epsilon\) for each pair, same shape as
r_ij(after broadcast).sum (bool, optional) – If True (default), return the sum over pairs; otherwise return per-pair energies.
- Returns:
Scalar total energy if
sumis True, else shape (P,) per-pair energies.- Return type:
- class torchff.vdw.Vdw(function: Literal['LennardJones', 'AmoebaVdw147'] = 'LennardJones', cutoff: float | None = None, use_customized_ops: bool = False, use_type_pairs: bool = False, sum_output: bool = True, cuda_graph_compat: bool = True)[source]#
Bases:
ModuleVan der Waals pair energy module (Lennard-Jones 12-6 or AMOEBA buffered 14-7).
Dispatches to
compute_lennard_jones_energy()/compute_vdw_14_7_energy()whenuse_customized_opsis True; otherwise uses minimum-image displacements viatorchff.pbc.PBCand the reference formulascompute_lennard_jones_energy_ref()/compute_vdw_14_7_energy_ref().- __init__(function: Literal['LennardJones', 'AmoebaVdw147'] = 'LennardJones', cutoff: float | None = None, use_customized_ops: bool = False, use_type_pairs: bool = False, sum_output: bool = True, cuda_graph_compat: bool = True)[source]#
- Parameters:
function ({'LennardJones', 'AmoebaVdw147'}, optional) – Potential form: standard LJ 12-6 or AMOEBA buffered 14-7.
cutoff (float, optional) – Stored on the module; the active cutoff is the
cutoffargument toforward().use_customized_ops (bool, optional) – If True, use custom CUDA/C++ kernels; otherwise use the PyTorch reference path.
use_type_pairs (bool, optional) – If True,
sigmaandepsilonare indexed byatom_typesfor each pair (shape(n_types, n_types)).sum_output (bool, optional) – If True (default), return a scalar sum over pairs. Must be True when
use_customized_opsis True because the custom kernels only return total energy. Whenuse_customized_opsis False, if False return per-pair energies of shape(P,).cuda_graph_compat (bool, optional) – If True (default), apply the cutoff with
torch.where()so tensor shapes are stable; if False, distances are filtered with boolean indexing before the energy expression.
- forward(coords: Tensor, pairs: Tensor, box: Tensor, sigma: Tensor, epsilon: Tensor, cutoff: float, atom_types: Tensor | None = None)[source]#
Compute vdW energy for the configured potential.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
pairs (torch.Tensor) – Shape (P, 2), pair indices (i, j).
box (torch.Tensor) – Periodic box, same convention as
torchff.pbc.PBC.sigma (torch.Tensor) – Per-pair
(P,)or type table(T, T)whenuse_type_pairsis True.epsilon (torch.Tensor) – Same layout as
sigma.cutoff (float) – Pair distance cutoff.
atom_types (torch.Tensor, optional) – Shape (N,), integer atom types; required when
use_type_pairsis True.
- Returns:
If
use_customized_opsis True, scalar total energy from the custom op. Otherwise per-pair energies of shape (P,), or a scalar ifsum_outputis True.- Return type:
Dispersion#
Tang–Tonnies damped C6 dispersion: U = -C6 f_6(br) / r^6.
- torchff.dispersion.compute_tang_tonnies_dispersion_energy(coords: Tensor, pairs: Tensor, box: Tensor, c6: Tensor, b: Tensor, cutoff: float, atom_types: Tensor | None = None) Tensor[source]#
Total Tang–Tonnies C6 dispersion energy (custom CUDA).
- torchff.dispersion.compute_tang_tonnies_dispersion_energy_ref(r_ij: Tensor, c6_ij: Tensor, b_ij: Tensor, sum: bool = True) Tensor[source]#
Reference: U_ij = -C6_ij f_6(b_ij r_ij) / r_ij^6.
- Parameters:
r_ij (torch.Tensor) – Pair distances, shape (P,) or broadcastable.
c6_ij (torch.Tensor) – Per-pair C6 and Tang–Tonnies inverse length b, same shape as
r_ijafter broadcast.b_ij (torch.Tensor) – Per-pair C6 and Tang–Tonnies inverse length b, same shape as
r_ijafter broadcast.sum (bool) – If True, return scalar sum over pairs; else per-pair energies (P,).
- class torchff.dispersion.Dispersion(cutoff: float | None = None, use_customized_ops: bool = False, use_type_pairs: bool = False, sum_output: bool = False, cuda_graph_compat: bool = True)[source]#
Bases:
ModuleTang–Tonnies C6 dispersion pair energy.
Dispatches to
compute_tang_tonnies_dispersion_energy()whenuse_customized_opsis True; otherwise PBC minimum-image distances andcompute_tang_tonnies_dispersion_energy_ref().
Slater#
Slater-type pair energies: A * P(B r) * exp(-B r) with P = (B r)^2 / 3 + B r + 1.
- torchff.slater.compute_slater_energy(coords: Tensor, pairs: Tensor, box: Tensor, A: Tensor, B: Tensor, cutoff: float, atom_types: Tensor | None = None) Tensor[source]#
Compute Slater pair energies via custom CUDA/C++ ops.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
pairs (torch.Tensor) – Shape (P, 2), integer indices (i, j) of interacting pairs.
box (torch.Tensor) – Shape (3, 3) or broadcastable, periodic box (same convention as
torchff.pbc).A (torch.Tensor) – Per-pair or type-pair amplitude \(A_{ij}\) (see
atom_types).B (torch.Tensor) – Per-pair or type-pair inverse length \(B_{ij}\) (see
atom_types).cutoff (float) – Distance cutoff; interactions beyond cutoff are excluded by the kernel.
atom_types (torch.Tensor, optional) – If provided, used by the backend for type-based indexing together with
AandB.
- Returns:
Scalar total Slater energy.
- Return type:
- torchff.slater.compute_slater_energy_ref(r_ij: Tensor, A_ij: Tensor, B_ij: Tensor, sum: bool = True)[source]#
Reference Slater pair energy in PyTorch.
With \(x = B_{ij} r_{ij}\), \(P = x^2/3 + x + 1\), \(E_{ij} = A_{ij} P \exp(-x)\).
- Parameters:
r_ij (torch.Tensor) – Pair distances, shape (P,) or broadcastable.
A_ij (torch.Tensor) – \(A_{ij}\) for each pair, same shape as
r_ij(after broadcast).B_ij (torch.Tensor) – \(B_{ij}\) for each pair, same shape as
r_ij(after broadcast).sum (bool, optional) – If True (default), return the sum over pairs; otherwise return per-pair energies.
- Returns:
Scalar total energy if
sumis True, else shape (P,) per-pair energies.- Return type:
- class torchff.slater.Slater(cutoff: float | None = None, use_customized_ops: bool = False, use_type_pairs: bool = False, sum_output: bool = False, cuda_graph_compat: bool = True)[source]#
Bases:
ModuleSlater pair energy module.
Dispatches to
compute_slater_energy()whenuse_customized_opsis True; otherwise uses minimum-image displacements viatorchff.pbc.PBCandcompute_slater_energy_ref().- __init__(cutoff: float | None = None, use_customized_ops: bool = False, use_type_pairs: bool = False, sum_output: bool = False, cuda_graph_compat: bool = True)[source]#
- Parameters:
cutoff (float, optional) – Stored on the module; the active cutoff is the
cutoffargument toforward().use_customized_ops (bool, optional) – If True, use custom CUDA/C++ kernels; otherwise use the PyTorch reference path.
use_type_pairs (bool, optional) – If True,
AandBare indexed byatom_typesfor each pair (shape(n_types, n_types)).sum_output (bool, optional) – If True, return a scalar sum over pairs. Requires
use_customized_opsFalse.cuda_graph_compat (bool, optional) – If True (default), apply the cutoff with
torch.where()so tensor shapes are stable; if False, distances are filtered with boolean indexing before the energy expression.
- forward(coords: Tensor, pairs: Tensor, box: Tensor, A: Tensor, B: Tensor, cutoff: float, atom_types: Tensor | None = None)[source]#
Compute Slater energy.
- Parameters:
coords (torch.Tensor) – Shape (N, 3), atom coordinates.
pairs (torch.Tensor) – Shape (P, 2), pair indices (i, j).
box (torch.Tensor) – Periodic box, same convention as
torchff.pbc.PBC.A (torch.Tensor) – Per-pair
(P,)or type table(T, T)whenuse_type_pairsis True.B (torch.Tensor) – Same layout as
A.cutoff (float) – Pair distance cutoff.
atom_types (torch.Tensor, optional) – Shape (N,), integer atom types; required when
use_type_pairsis True.
- Returns:
If
use_customized_opsis True, scalar total energy from the custom op. Otherwise per-pair energies of shape (P,), or a scalar ifsum_outputis True.- Return type:
Neighbor List#
- class torchff.nblist.NeighborList(natoms: int, exclusions: Tensor | None = None, include_self: bool = False, use_customized_ops: bool = True, algorithm: str = 'nsquared')[source]#
Bases:
Module- __init__(natoms: int, exclusions: Tensor | None = None, include_self: bool = False, use_customized_ops: bool = True, algorithm: str = 'nsquared')[source]#
- classmethod convert_pairs_coo_to_csr(pairs_coo: Tensor, num_atoms: int | None = None) Tuple[Tensor, Tensor][source]#
Convert a COO pair list into CSR format for a symmetric matrix.
Given a list of index pairs in COO (coordinate) format representing non-zero entries of a symmetric matrix, this method produces the equivalent CSR (Compressed Sparse Row) representation. Because the matrix is symmetric, both directions
(i, j)and(j, i)are guaranteed in the output regardless of whether the input contains one or both directions.The conversion proceeds in four steps:
Symmetrize – concatenate
pairs_coowith its column-flipped copy so that both(i, j)and(j, i)are present.Encode and deduplicate – map each pair to a scalar key
row * num_atoms + coland calltorch.unique()on the 1-D key tensor. This is more efficient thantorch.unique(dim=0)and produces keys in sorted (row-major) order.Decode – recover row and column indices from the unique keys via integer division and modulo.
Build
row_ptr– count neighbors per row withtorch.Tensor.scatter_add_(), then take the cumulative sum to form the CSR indptr array.
- Parameters:
pairs_coo (torch.Tensor) – Integer tensor of shape
(N, 2)where each row[i, j]denotes a non-zero entry. The input may contain only one direction of a symmetric pair, both directions, or a mix.num_atoms (int, optional) – Total number of atoms (i.e. the matrix dimension). When
None(default), it is inferred aspairs_coo.max() + 1. Specify this explicitly when isolated atoms with no pairs exist so thatrow_ptrhas the correct length.
- Returns:
row_ptr (torch.Tensor) – Shape
(num_atoms + 1,), dtypetorch.long.row_ptr[i]is the start index incol_indicesfor row i;row_ptr[i+1] - row_ptr[i]gives the number of neighbors of atom i.col_indices (torch.Tensor) – Shape
(M,), dtypetorch.long. Column indices of all non-zero entries, sorted first by row then by column within each row.Mis the total number of non-zero entries after symmetrization and deduplication.
- torchff.nblist.build_neighbor_list_nsquared(coords: Tensor, box: Tensor, cutoff: float, max_npairs: int = -1, padding: bool = False, excl_row_ptr: Tensor | None = None, excl_col_indices: Tensor | None = None, include_self: bool = False, out: Tensor | None = None)[source]#
Build a neighbor list using the O(N^2) algorithm.
- Parameters:
coords (torch.Tensor) – Atom positions, shape
(natoms, 3).box (torch.Tensor) – Periodic box vectors, shape
(3, 3).cutoff (float) – Distance cutoff.
max_npairs (int, optional) – Pre-allocated capacity for pairs.
-1allocates for the worst casenatoms*(natoms-1)/2.padding (bool, optional) – If
True, return the full (possibly padded) pairs tensor instead of trimming to the actual count.excl_row_ptr (torch.Tensor, optional) – CSR row-pointer tensor of shape
(natoms + 1,)for the exclusion list. Must betorch.long.excl_col_indices (torch.Tensor, optional) – CSR column-indices tensor for the exclusion list. Must be
torch.long. Pairs(i, j)present in this sparse structure are skipped during neighbor-list construction. UseNeighborList.convert_pairs_coo_to_csr()to convert COO-format exclusion pairs into the required CSR tensors.include_self (bool, optional) – If
True, include self-pairs(i, i).out (torch.Tensor, optional) – Pre-allocated output tensor of shape
(max_npairs, 2)(dtypetorch.long). When provided, the kernel writes into this tensor andmax_npairsis ignored.
- Returns:
pairs (torch.Tensor) – Neighbor pairs, shape
(npairs, 2).npairs (torch.Tensor) – Scalar tensor with the number of pairs found.
- torchff.nblist.build_neighbor_list_cell_list(coords: Tensor, box: Tensor, cutoff: float, max_npairs: int = -1, cell_size: float = 0.4, padding: bool = False, shared: bool = False)[source]#
- torchff.nblist.build_cluster_pairs(coords: Tensor, box: Tensor, cutoff: float, exclusions: Tensor | None = None, cell_size: float = 0.4, max_num_interacting_clusters: int = -1) Tuple[Tensor, Tensor, Tensor, Tensor][source]#
Nonbonded#
- torchff.nonbonded.compute_nonbonded_energy_from_atom_pairs(coords: Tensor, pairs: Tensor, box: Tensor, sigma: Tensor, epsilon: Tensor, charges: Tensor, coul_constant, cutoff: float, do_shift: bool = True) Tensor[source]#
Compute fused Coulomb + Lennard-Jones nonbonded energies using custom CUDA/C++ ops.
- torchff.nonbonded.compute_nonbonded_energy_from_atom_pairs_ref(coords: Tensor, pairs: Tensor, box: Tensor, sigma: Tensor, epsilon: Tensor, charges: Tensor, coul_constant: float, cutoff: float, do_shift: bool = True) Tensor[source]#
Reference fused Coulomb + Lennard-Jones implementation using native PyTorch ops.
- torchff.nonbonded.compute_nonbonded_forces_from_atom_pairs(coords: Tensor, pairs: Tensor, box: Tensor, sigma: Tensor, epsilon: Tensor, charges: Tensor, coul_constant: float, cutoff: float, forces: Tensor) Tensor[source]#
Compute fused Coulomb + Lennard-Jones nonbonded forces in-place using custom CUDA/C++ ops.
Multipoles#
- torchff.multipoles.computeInteractionTensor(drVec: Tensor, dampFactors: Tensor | None = None, drInv: Tensor | None = None, rank: int = 2)[source]#
- class torchff.multipoles.MultipolePacker(rank: int = 2)[source]#
Bases:
ModuleConvert monopole, dipole, and quadrupole to (N, 10) polytensor with quadrupole entries scaled so symmetry-equivalent operations are avoided.
- forward(mono: Tensor, dipo: Tensor | None = None, quad: Tensor | None = None) Tensor[source]#
- Parameters:
mono (torch.Tensor) – Monopoles, shape (N,).
dipo (torch.Tensor) – Dipoles, shape (N, 3).
quad (torch.Tensor) – Quadrupoles, shape (N, 3, 3).
- Returns:
Polytensor [q, ux, uy, uz, Qxx, Qxy, Qxz, Qyy, Qyz, Qzz], shape (N, 10).
- Return type:
- class torchff.multipoles.MultipolarInteraction(rank: int, cutoff: float, ewald_alpha: float = -1.0, prefactor: float = 1.0, return_fields: bool = False, use_customized_ops: bool = True, cuda_graph_compat: bool = True)[source]#
Bases:
Module
- class torchff.multipoles.AxisTypes(*values)[source]#
Bases:
IntEnum- ZThenX = 0#
- Bisector = 1#
- ZBisect = 2#
- ThreeFold = 3#
- ZOnly = 4#
- NoAxisType = 5#
- LastAxisTypeIndex = 6#
- class torchff.multipoles.MultipolarRotation(use_customized_ops: bool = False)[source]#
Bases:
ModuleMultipole local frame: build rotation matrices and rotate dipole / quadrupole tensors.
Vectorized local-to-global rotation matrices for multipole sites use the same batched construction as the historical TorchFF Python path: differences
positions[neighbor] - positions[site]with no periodic minimum-image wrapping. Intramolecular (or otherwise local) geometry should already lie in one periodic image so that these raw vectors match the intended local frame.- Parameters:
use_customized_ops (bool) – If
True, use_compute_rotation_matrices_torchff()(CUDA custom op); ifFalse, use_compute_rotation_matrices_python()(Python reference). Dipole and quadrupole rotation userotateDipoles()androtateQuadrupoles()in either case.
Notes
Rotation matrices have shape
(N, 3, 3); rows are local X, Y, Z in global coordinates (hstackofxVec,yVec,zVec). Inputsz_atoms,x_atoms,y_atoms, andaxis_typeshave shape(N,)(seeAxisTypes).- classmethod compute_matrices(coords: Tensor, z_atoms: Tensor, x_atoms: Tensor, y_atoms: Tensor, axis_types: Tensor, *, use_customized_ops: bool = False) Tensor[source]#
Shape
(N, 3, 3)rotation matrices (rows = local X, Y, Z in global coordinates).
- classmethod rotate_dipoles(matrices: Tensor, dipoles: Tensor) Tensor[source]#
Rotate dipoles of shape
(N, 3); returns(N, 1, 3)(same asrotateDipoles()).
- classmethod rotate_quadrupoles(matrices: Tensor, quadrupoles: Tensor) Tensor[source]#
Rotate quadrupoles of shape
(N, 3, 3).
- forward(coords: Tensor, z_atoms: Tensor, x_atoms: Tensor, y_atoms: Tensor, axis_types: Tensor, dipoles: Tensor, quadrupoles: Tensor | None = None) Tensor | Tuple[Tensor, Tensor][source]#
Apply local-to-global rotation to dipoles and optionally quadrupoles.
- Returns:
Rotated dipoles
(N, 3). Ifquadrupolesis given, returns(dipoles_rot, quadrupoles_rot).- Return type:
- torchff.multipoles.scaleMultipoles(mPoles: Tensor, monoScales: Tensor, dipoScales: Tensor, quadScales: Tensor)[source]#
- torchff.multipoles.rotateMultipoles(mono: Tensor, dipo: Tensor, quad: Tensor, rotMatrix: Tensor)[source]#
Rotate multipoles
- Parameters:
mono (torch.Tensor) – Monopoles, shape (N,)
dipo (torch.Tensor) – Dipoles, shape (N, 3)
quad (torch.Tensor) – Quadrupoles, shape (N, 3, 3)
- Returns:
mPoles – Multipoles [q, ux, uy, uz, Qxx, Qxy, Qxz, Qyy, Qyz, Qzz], shape (N, 10)
- Return type:
- torchff.multipoles.convertMultipolesToPolytensor(mono: Tensor, dipo: Tensor, quad: Tensor)[source]#
Takes already-rotated multipoles and flattens to (N, 10) polytensor with quadrupole entries appropriately scaled so that symmetry-equivalent operations are avoided.
- Parameters:
mono (torch.Tensor) – Monopoles, shape (N,)
dipo (torch.Tensor) – Dipoles, shape (N, 3)
quad (torch.Tensor) – Quadrupoles, shape (N, 3, 3)
- Returns:
mPoles – Multipoles [q, ux, uy, uz, Qxx, Qxy, Qxz, Qyy, Qyz, Qzz], shape (N, 10)
- Return type:
- torchff.multipoles.computeCartesianQuadrupoles(quad_s: Tensor)[source]#
Compute cartesian quadrupoles from spheric-harmonics quadrupoles
- Parameters:
quad_s (torch.Tensor) – Quadrupoles in spherical harmonics form (Q20, Q21c, Q21s, Q22c, Q22s), shape (N, 5).
- Returns:
quad – Quadrupoles in cartesian form, shape N x 3 x 3
- Return type:
- torchff.multipoles.computeSphericalQuadrupoles(quad_c: Tensor)[source]#
Compute cartesian quadrupoles from spheric-harmonics quadrupoles
- Parameters:
quad_c (torch.Tensor) – Quadrupoles in cartesian form (Qxx, Qxy, Qxz, Qyy, Qyz, Qzz), shape (N, 6).
- Returns:
quad_s – Quadrupoles in spherical harmonics form, shape (N, 5)
- Return type:
Ewald Summation#
- torchff.ewald.ewald_long_range(coords: Tensor, box: Tensor, q: Tensor, p: Tensor | None, t: Tensor | None, kmax: int, alpha: float)[source]#
- torchff.ewald.ewald_long_range_with_fields(coords: Tensor, box: Tensor, q: Tensor, p: Tensor | None, t: Tensor | None, kmax: int, alpha: float)[source]#
Particle Mesh Ewald#
PME Reference#
B-Splines#
B-spline utilities for PME, matching the logic in csrc/pme/bsplines.cuh.
- torchff.bsplines.get_bspline_coeff_order6(i: int) float[source]#
B-spline coefficient for order-6 PME at index i.
Matches
get_bspline_coeff_order6()incsrc/pme/bsplines.cuh.
- torchff.bsplines.get_bspline_modulus(k: Tensor, K: int, order: int = 6) Tensor[source]#
B-spline modulus for PME reciprocal-space scaling.
For each wave vector index k, computes sum_m b(m) * cos(2*pi*m*k/K) where m runs over the order-6 B-spline range and b(m) are the B-spline coefficients. Matches
get_bspline_modulus_device()incsrc/pme/bsplines.cuh.- Parameters:
k (torch.Tensor) – 1D tensor of integer indices (e.g. 0 .. K-1 or 0 .. K//2).
K (int) – Grid size along this dimension.
order (int) – B-spline order; only 6 is implemented.
- Returns:
Modulus values, same shape and device as k, dtype float.
- Return type:
- torchff.bsplines.compute_bspline_moduli_1d(K: int, dtype: dtype = torch.float32, device: device = None) Tensor[source]#
Precompute the 1D B-spline modulus array for grid size K.
Returns moduli for indices 0, 1, …, K-1. For use as xmoduli or ymoduli in PME (length K).
- Parameters:
K (int) – Grid size along this dimension.
dtype (torch.dtype) – Output tensor dtype.
device (torch.device, optional) – Output device; default is CPU.
- Returns:
Shape (K,), dtype and device as requested.
- Return type:
- torchff.bsplines.compute_bspline_moduli_z(K3: int, dtype: dtype = torch.float32, device: device = None) Tensor[source]#
Precompute the z-dimension B-spline modulus array for real FFT.
Returns moduli for z indices 0, 1, …, K3//2 (length K3//2 + 1), as used in the conjugate-even (rfft) reciprocal grid.
- Parameters:
K3 (int) – Grid size along the z dimension.
dtype (torch.dtype) – Output tensor dtype.
device (torch.device, optional) – Output device; default is CPU.
- Returns:
Shape (K3//2 + 1,), dtype and device as requested.
- Return type:
Periodic Boundary Conditions#
- class torchff.pbc.PBC(box: Tensor | None = None, box_inv: Tensor | None = None)[source]#
Bases:
ModuleApply periodic boundary conditions to vectors.
Supports both fixed box (provided at construction) and variable box (passed to
forward()). If neither box nor box_inv is available, returns the input vectors unchanged.- forward(dr_vecs: Tensor, box: Tensor | None = None, box_inv: Tensor | None = None) Tensor[source]#
Apply PBC to displacement vectors.
- Parameters:
dr_vecs (torch.Tensor) – Real space vectors in Cartesian, shape (N, 3).
box (torch.Tensor, optional) – Simulation box with axes in rows, shape (3, 3). Overrides the buffer if provided.
box_inv (torch.Tensor, optional) – Inverse of the box matrix, shape (3, 3). Overrides the buffer if provided.
- Returns:
PBC-folded vectors in Cartesian, shape (N, 3).
- Return type:
- torchff.pbc.applyPBC(dr_vecs: Tensor, box: Tensor | None = None, box_inv: Tensor | None = None) Tensor[source]#
Apply periodic boundary conditions to a set of vectors.
Functional interface to
PBC. Prefer usingPBCwhen integrating withnn.Module-based models.- Parameters:
dr_vecs (torch.Tensor) – Real space vectors in Cartesian, shape (N, 3).
box (torch.Tensor, optional) – Simulation box with axes in rows, shape (3, 3).
box_inv (torch.Tensor, optional) – Inverse of the box matrix, shape (3, 3).
- Returns:
PBC-folded vectors in Cartesian, shape (N, 3).
- Return type: