Crate potentials

Crate potentials 

Source
Expand description

§Potentials

A high-performance, pure Rust library for classical molecular dynamics potential energy functions.

Designed for force field engines, MD simulators, and scientific computing in embedded environments.

FeaturesInstallationUsageModulesPerformanceLicense


§Features

  • 🚀 High Performance

    • Zero Heap Allocation: All computations happen on the stack
    • Branchless Kernels: Mask-based conditionals for vectorization-friendly code
    • Optimized Implementations: Shared subexpressions, Horner’s method, Chebyshev recursion
    • ~650M interactions/sec: LJ energy+force throughput on modern CPUs
  • 🎯 Comprehensive Coverage

    • 7 Pair Potentials: LJ, Mie, Buckingham, Coulomb, Yukawa, Gaussian, Soft sphere
    • 6 Bond Potentials: Harmonic, GROMOS96, Morse, FENE, Cubic, Quartic
    • 5 Angle Potentials: Harmonic, Cosine, Linear, Urey-Bradley, Cross
    • 4 Torsion Potentials: Periodic cosine, OPLS, Ryckaert-Bellemans, Harmonic
    • 3 Improper Potentials: Harmonic, Cosine, Distance-based
    • 2 H-Bond Potentials: DREIDING 12-10, LJ-Cosine
    • 6 Meta Wrappers: Cutoff, Shift, Switch, Softcore, Scaled, Sum
  • 🔬 Physically Accurate

    • Analytically correct energy and force expressions
    • Numerical derivative validation for every potential
    • Proper handling of edge cases and singularities
    • CODATA 2018 physical constants
  • 🛰️ Embeddable & Portable

    • no_std Support: Works in embedded devices, kernels, or WASM
    • Pure Rust: No C/C++ dependencies, simple build chain
    • Generic Precision: Supports both f32 and f64

§Installation

Add this to your Cargo.toml:

[dependencies]
potentials = "0.1.0"

To use in a no_std environment:

[dependencies]
potentials = { version = "0.1.0", default-features = false, features = ["libm"] }

§Usage

§Quick Start

use potentials::{pair::Lj, base::Potential2};

// Create Lennard-Jones potential (ε=1.0 kcal/mol, σ=3.4 Å)
let lj = Lj::<f64>::new(1.0, 3.4);

// Compute energy and force factor at r=4.0 Å
let r_sq = 16.0;  // r² = 4²
let (energy, force_factor) = lj.energy_force(r_sq);

// Force vector: F_ij = force_factor * r_ij_vec
println!("Energy: {:.6} kcal/mol", energy);
println!("Force factor: {:.6}", force_factor);

§Adding Cutoffs and Modifications

use potentials::{pair::Lj, meta::Switch, base::Potential2};

let lj = Lj::<f64>::new(1.0, 3.4);

// Smooth switching function
let lj_switched: Switch<_, f64> = Switch::new(lj, 9.0, 12.0);

// Energy smoothly goes to zero between r=9 and r=12 Å
let energy = lj_switched.energy(100.0);  // r² = 100

§Combining Potentials

use potentials::{pair::{Lj, Coul}, meta::Sum, base::Potential2};

let lj = Lj::<f64>::new(0.1, 3.0);
let coul = Coul::<f64>::new(-332.0636);  // Na⁺-Cl⁻ attraction

// LJ + Coulomb combined potential
let combined: Sum<_, _, f64> = Sum::new(lj, coul);

let (energy, force_factor) = combined.energy_force(25.0);  // r² = 25

§Angle and Torsion Potentials

use potentials::{angle::Cos, torsion::Opls, base::{Potential3, Potential4}};

// Water H-O-H angle (cosine form, faster than harmonic)
let angle = Cos::<f64>::from_degrees(100.0, 104.5);
let (e_angle, d_angle) = angle.energy_derivative(1.0, 1.0, 0.25);

// Alkane torsion (OPLS Fourier series)
let torsion = Opls::<f64>::new(0.7, 0.1, 1.5, 0.0);
let (e_tors, d_tors) = torsion.energy_derivative(0.5, 0.866);

§Modules

§Architecture

The library is organized around three core traits:

TraitBodiesInputOutputUse Case
Potential22(V, S)Pairs, Bonds
Potential33r_ij², r_jk², cos θ(V, dV/d(cos θ))Angles
Potential44cos φ, sin φ(V, dV/dφ)Torsions

§Module Reference

ModulePotentialsDescription
pairLj, Mie, Buck, Coul, Yukawa, Gauss, SoftNon-bonded pair interactions
bondHarm, G96, Morse, Fene, Cubic, QuartCovalent bond stretching
angleHarm, Cos, Linear, Urey, CrossValence angle bending
torsionCos, Opls, Rb, HarmProper dihedral angles
impHarm, Cos, DistImproper torsions
hbondDreid, LjCosHydrogen bonding
metaCutoff, Shift, Switch, Softcore, Scaled, SumPotential modifiers
constsPhysical constants (CODATA 2018)

§Force Field Compatibility

Force FieldPairBondAngleTorsionImproper
AMBER
CHARMM
OPLS
GROMOS
DREIDING

§Performance

Benchmarked on an Intel Core i7-13620H Laptop CPU (10 cores).

§Pair Potentials (energy + force)

PotentialThroughputTime per eval
Lennard-Jones200 M/s5.0 ns
Coulomb69 M/s14.5 ns
Gaussian57 M/s17.6 ns
Soft sphere53 M/s19.0 ns
Yukawa41 M/s24.4 ns
Buckingham36 M/s27.4 ns
Mie (n-m)28 M/s36.3 ns

§Angle Potentials

PotentialThroughputNotes
Cosine367 M/sNo acos needed
Harmonic54 M/sRequires acos

§Meta Wrapper Overhead

ConfigurationThroughputOverhead
LJ bare200 M/sbaseline
LJ + softcore113 M/s1.8×
LJ + cutoff107 M/s1.9×
LJ + shift100 M/s2.0×
LJ + switch76 M/s2.6×

§Throughput Scaling

InteractionsThroughputNotes
1,000660 M/sCache-hot
10,000648 M/sL2 cache
100,000649 M/sL3 cache

§License

This project is licensed under the MIT License - see the LICENSE file for details.


Made with ❤️ for the scientific computing community

Modules§

angle
Angle Bending Potentials
base
Physics Interface Traits
bond
Bonded Stretch Potentials
consts
Physical Constants
hbond
Hydrogen Bond Potentials
imp
Improper Torsion Potentials
math
Math Abstraction Layer
meta
Meta Potentials (Wrappers and Modifiers)
pair
Non-Bonded Pair Potentials
torsion
Dihedral Torsion Potentials