Skip to main content

sci_form/ir/
vibrations.rs

1//! Vibrational analysis: normal modes, frequencies, IR intensities, and spectrum.
2//!
3//! Given a numerical Hessian, computes mass-weighted normal modes,
4//! vibrational frequencies (cm⁻¹), IR intensities from dipole derivatives,
5//! and generates a Lorentzian-broadened IR spectrum.
6
7use nalgebra::DMatrix;
8use serde::{Deserialize, Serialize};
9
10use super::hessian::{compute_numerical_hessian, HessianMethod};
11
12/// Atomic masses in amu for elements 1–86 (H through Rn).
13fn atomic_mass(z: u8) -> f64 {
14    match z {
15        1 => 1.00794,
16        2 => 4.00260,
17        5 => 10.811,
18        6 => 12.0107,
19        7 => 14.0067,
20        8 => 15.9994,
21        9 => 18.9984,
22        14 => 28.0855,
23        15 => 30.9738,
24        16 => 32.065,
25        17 => 35.453,
26        22 => 47.867,
27        24 => 51.9961,
28        25 => 54.9380,
29        26 => 55.845,
30        27 => 58.9332,
31        28 => 58.6934,
32        29 => 63.546,
33        30 => 65.38,
34        35 => 79.904,
35        44 => 101.07,
36        46 => 106.42,
37        47 => 107.868,
38        53 => 126.904,
39        78 => 195.084,
40        79 => 196.967,
41        _ => z as f64 * 1.5, // rough fallback
42    }
43}
44
45/// A single vibrational normal mode.
46#[derive(Debug, Clone, Serialize, Deserialize)]
47pub struct VibrationalMode {
48    /// Frequency in cm⁻¹ (negative = imaginary).
49    pub frequency_cm1: f64,
50    /// IR intensity in km/mol.
51    pub ir_intensity: f64,
52    /// Displacement vector (3N elements, mass-weighted normal coordinate).
53    pub displacement: Vec<f64>,
54    /// Whether this is a real mode (true) or translation/rotation (false).
55    pub is_real: bool,
56}
57
58/// Complete vibrational analysis result.
59#[derive(Debug, Clone, Serialize, Deserialize)]
60pub struct VibrationalAnalysis {
61    /// Number of atoms.
62    pub n_atoms: usize,
63    /// All vibrational modes sorted by frequency.
64    pub modes: Vec<VibrationalMode>,
65    /// Number of real vibrational modes (excluding translation/rotation).
66    pub n_real_modes: usize,
67    /// Zero-point vibrational energy in eV.
68    pub zpve_ev: f64,
69    /// Semiempirical method used for Hessian.
70    pub method: String,
71    /// Notes and caveats.
72    pub notes: Vec<String>,
73}
74
75/// A peak in the IR spectrum.
76#[derive(Debug, Clone, Serialize, Deserialize)]
77pub struct IrPeak {
78    /// Frequency in cm⁻¹.
79    pub frequency_cm1: f64,
80    /// IR intensity in km/mol.
81    pub ir_intensity: f64,
82    /// Mode index.
83    pub mode_index: usize,
84}
85
86/// Broadened IR spectrum.
87#[derive(Debug, Clone, Serialize, Deserialize)]
88pub struct IrSpectrum {
89    /// Wavenumber axis (cm⁻¹).
90    pub wavenumbers: Vec<f64>,
91    /// Transmittance/absorbance axis.
92    pub intensities: Vec<f64>,
93    /// Discrete peaks from vibrational analysis.
94    pub peaks: Vec<IrPeak>,
95    /// Broadening width in cm⁻¹.
96    pub gamma: f64,
97    /// Notes.
98    pub notes: Vec<String>,
99}
100
101// Constants for unit conversion
102// 1 eV = 8065.544 cm⁻¹
103const EV_TO_CM1: f64 = 8065.544;
104// 1 amu = 1.66054e-27 kg, 1 eV = 1.60218e-19 J, 1 Å = 1e-10 m
105// Hessian in eV/Ų, mass in amu → frequency in cm⁻¹:
106// ω² = H/(m) in eV/(amu·Å²) → convert to s⁻²:
107// factor = eV_to_J / (amu_to_kg * Å_to_m²) = 1.60218e-19 / (1.66054e-27 * 1e-20) = 9.6485e27
108// ν_cm1 = sqrt(factor * eigenvalue) / (2π * c) where c = 2.998e10 cm/s
109const HESSIAN_TO_FREQUENCY_FACTOR: f64 = 9.6485e27; // eV/(amu·Å²) → s⁻²
110const INV_2PI_C: f64 = 1.0 / (2.0 * std::f64::consts::PI * 2.99792458e10); // 1/(2πc) in s/cm
111
112/// Compute the molecular dipole for use in IR intensity calculation.
113fn compute_dipole_vector(
114    elements: &[u8],
115    positions: &[[f64; 3]],
116    method: HessianMethod,
117) -> Result<[f64; 3], String> {
118    match method {
119        HessianMethod::Eht => {
120            let eht = crate::eht::solve_eht(elements, positions, None)?;
121            let dipole = crate::dipole::compute_dipole_from_eht(
122                elements,
123                positions,
124                &eht.coefficients,
125                eht.n_electrons,
126            );
127            Ok(dipole.vector)
128        }
129        HessianMethod::Pm3 => {
130            let pm3 = crate::pm3::solve_pm3(elements, positions)?;
131            // Use Mulliken charges for dipole
132            let dipole = crate::dipole::compute_dipole(&pm3.mulliken_charges, positions);
133            Ok(dipole.vector)
134        }
135        HessianMethod::Xtb => {
136            let xtb = crate::xtb::solve_xtb(elements, positions)?;
137            let dipole = crate::dipole::compute_dipole(&xtb.mulliken_charges, positions);
138            Ok(dipole.vector)
139        }
140    }
141}
142
143/// Compute IR intensities from numerical dipole derivatives along normal modes.
144///
145/// I_k ∝ |dμ/dQ_k|² where Q_k is the k-th normal coordinate.
146fn compute_ir_intensities(
147    elements: &[u8],
148    positions: &[[f64; 3]],
149    normal_modes: &DMatrix<f64>,
150    method: HessianMethod,
151    delta: f64,
152) -> Result<Vec<f64>, String> {
153    let n_atoms = elements.len();
154    let n_modes = normal_modes.ncols();
155    let masses: Vec<f64> = elements.iter().map(|&z| atomic_mass(z)).collect();
156
157    let mut intensities = Vec::with_capacity(n_modes);
158
159    for k in 0..n_modes {
160        // Displace along normal mode k
161        let mut pos_plus: Vec<[f64; 3]> = positions.to_vec();
162        let mut pos_minus: Vec<[f64; 3]> = positions.to_vec();
163
164        for i in 0..n_atoms {
165            let sqrt_m = masses[i].sqrt();
166            for c in 0..3 {
167                let disp = normal_modes[(3 * i + c, k)] * delta / sqrt_m;
168                pos_plus[i][c] += disp;
169                pos_minus[i][c] -= disp;
170            }
171        }
172
173        let mu_plus = compute_dipole_vector(elements, &pos_plus, method)?;
174        let mu_minus = compute_dipole_vector(elements, &pos_minus, method)?;
175
176        // dμ/dQ ≈ (μ+ - μ-) / (2δ)
177        let dmu_dq: Vec<f64> = (0..3)
178            .map(|c| (mu_plus[c] - mu_minus[c]) / (2.0 * delta))
179            .collect();
180
181        // I ∝ |dμ/dQ|² — convert to km/mol with empirical scaling
182        let intensity = dmu_dq.iter().map(|x| x * x).sum::<f64>();
183        // Scale: Debye²/Ų·amu → km/mol (empirical: ~42.256)
184        intensities.push(intensity * 42.256);
185    }
186
187    Ok(intensities)
188}
189
190/// Lorentzian line shape: L(x) = (1/π) · γ / [(x - x₀)² + γ²]
191fn lorentzian(x: f64, x0: f64, gamma: f64) -> f64 {
192    let g = gamma.max(1e-6);
193    let dx = x - x0;
194    g / (std::f64::consts::PI * (dx * dx + g * g))
195}
196
197/// Perform a complete vibrational analysis from elements and positions.
198///
199/// Computes the numerical Hessian, mass-weighted eigenvalue problem,
200/// vibrational frequencies, and IR intensities.
201pub fn compute_vibrational_analysis(
202    elements: &[u8],
203    positions: &[[f64; 3]],
204    method: HessianMethod,
205    step_size: Option<f64>,
206) -> Result<VibrationalAnalysis, String> {
207    let n_atoms = elements.len();
208    let n3 = 3 * n_atoms;
209    let delta = step_size.unwrap_or(0.005);
210
211    if n_atoms < 2 {
212        return Err("Need at least 2 atoms for vibrational analysis".to_string());
213    }
214
215    // 1. Compute numerical Hessian
216    let hessian = compute_numerical_hessian(elements, positions, method, Some(delta))?;
217
218    // 2. Build mass-weighted Hessian: H'_{ij} = H_{ij} / sqrt(m_i * m_j)
219    let masses: Vec<f64> = elements.iter().map(|&z| atomic_mass(z)).collect();
220    let mut mw_hessian = DMatrix::zeros(n3, n3);
221    for i in 0..n3 {
222        let mi = masses[i / 3];
223        for j in 0..n3 {
224            let mj = masses[j / 3];
225            mw_hessian[(i, j)] = hessian[(i, j)] / (mi * mj).sqrt();
226        }
227    }
228
229    // 3. Diagonalize mass-weighted Hessian
230    let eigen = mw_hessian.symmetric_eigen();
231
232    // 4. Sort eigenvalues and eigenvectors by eigenvalue
233    let mut indices: Vec<usize> = (0..n3).collect();
234    indices.sort_by(|&a, &b| {
235        eigen.eigenvalues[a]
236            .partial_cmp(&eigen.eigenvalues[b])
237            .unwrap_or(std::cmp::Ordering::Equal)
238    });
239
240    // 5. Convert eigenvalues to frequencies (cm⁻¹)
241    // For linear molecules: 5 zero modes, for nonlinear: 6
242    let is_linear = n_atoms == 2; // simplified check
243    let _n_tr = if is_linear { 5 } else { 6 };
244    let freq_threshold = 50.0; // cm⁻¹ threshold for "real" modes
245
246    let mut sorted_eigenvalues = Vec::with_capacity(n3);
247    let mut sorted_modes = DMatrix::zeros(n3, n3);
248    for (new_idx, &old_idx) in indices.iter().enumerate() {
249        sorted_eigenvalues.push(eigen.eigenvalues[old_idx]);
250        for i in 0..n3 {
251            sorted_modes[(i, new_idx)] = eigen.eigenvectors[(i, old_idx)];
252        }
253    }
254
255    // Convert eigenvalues to frequencies
256    let frequencies: Vec<f64> = sorted_eigenvalues
257        .iter()
258        .map(|&ev| {
259            if ev >= 0.0 {
260                // ν = (1/2πc) * sqrt(eigenvalue * factor)
261                (ev * HESSIAN_TO_FREQUENCY_FACTOR).sqrt() * INV_2PI_C
262            } else {
263                // Imaginary frequency (negative eigenvalue)
264                -((-ev) * HESSIAN_TO_FREQUENCY_FACTOR).sqrt() * INV_2PI_C
265            }
266        })
267        .collect();
268
269    // 6. Compute IR intensities from dipole derivatives
270    let ir_intensities = compute_ir_intensities(
271        elements,
272        positions,
273        &sorted_modes,
274        method,
275        delta * 2.0, // slightly larger step for dipole derivatives
276    )?;
277
278    // 7. Build mode list
279    let mut modes = Vec::with_capacity(n3);
280    let mut n_real = 0;
281    let mut zpve = 0.0;
282
283    for k in 0..n3 {
284        let freq = frequencies[k];
285        let is_real = freq.abs() > freq_threshold;
286        if is_real && freq > 0.0 {
287            n_real += 1;
288            // ZPVE contribution: (1/2)hν in eV
289            zpve += 0.5 * freq / EV_TO_CM1;
290        }
291
292        let displacement: Vec<f64> = (0..n3).map(|i| sorted_modes[(i, k)]).collect();
293
294        modes.push(VibrationalMode {
295            frequency_cm1: freq,
296            ir_intensity: ir_intensities.get(k).copied().unwrap_or(0.0),
297            displacement,
298            is_real,
299        });
300    }
301
302    let method_name = match method {
303        HessianMethod::Eht => "EHT",
304        HessianMethod::Pm3 => "PM3",
305        HessianMethod::Xtb => "xTB",
306    };
307
308    Ok(VibrationalAnalysis {
309        n_atoms,
310        modes,
311        n_real_modes: n_real,
312        zpve_ev: zpve,
313        method: method_name.to_string(),
314        notes: vec![
315            format!(
316                "Numerical Hessian computed with {} using central finite differences (δ = {} Å).",
317                method_name, delta
318            ),
319            "IR intensities derived from numerical dipole derivatives along normal coordinates."
320                .to_string(),
321            "Frequencies below 50 cm⁻¹ are classified as translation/rotation modes.".to_string(),
322        ],
323    })
324}
325
326/// Generate a Lorentzian-broadened IR spectrum from vibrational analysis.
327///
328/// `gamma`: Lorentzian half-width in cm⁻¹ (typically 10–30).
329/// `wn_min`, `wn_max`: wavenumber range in cm⁻¹.
330/// `n_points`: number of grid points.
331pub fn compute_ir_spectrum(
332    analysis: &VibrationalAnalysis,
333    gamma: f64,
334    wn_min: f64,
335    wn_max: f64,
336    n_points: usize,
337) -> IrSpectrum {
338    let n_points = n_points.max(2);
339    let step = (wn_max - wn_min) / (n_points as f64 - 1.0);
340    let wavenumbers: Vec<f64> = (0..n_points).map(|i| wn_min + step * i as f64).collect();
341    let mut intensities = vec![0.0; n_points];
342
343    let mut peaks = Vec::new();
344
345    for (mode_idx, mode) in analysis.modes.iter().enumerate() {
346        if !mode.is_real || mode.frequency_cm1 <= 0.0 {
347            continue;
348        }
349
350        peaks.push(IrPeak {
351            frequency_cm1: mode.frequency_cm1,
352            ir_intensity: mode.ir_intensity,
353            mode_index: mode_idx,
354        });
355
356        for (idx, &wn) in wavenumbers.iter().enumerate() {
357            intensities[idx] += mode.ir_intensity * lorentzian(wn, mode.frequency_cm1, gamma);
358        }
359    }
360
361    peaks.sort_by(|a, b| {
362        b.ir_intensity
363            .partial_cmp(&a.ir_intensity)
364            .unwrap_or(std::cmp::Ordering::Equal)
365    });
366
367    IrSpectrum {
368        wavenumbers,
369        intensities,
370        peaks,
371        gamma,
372        notes: vec![
373            format!(
374                "IR spectrum generated with Lorentzian broadening (γ = {} cm⁻¹).",
375                gamma
376            ),
377            format!("Vibrational analysis method: {}.", analysis.method),
378        ],
379    }
380}
381
382#[cfg(test)]
383mod tests {
384    use super::*;
385
386    #[test]
387    fn test_lorentzian_peak_at_center() {
388        let val = lorentzian(100.0, 100.0, 10.0);
389        // At center: L = 1/(π·γ)
390        let expected = 1.0 / (std::f64::consts::PI * 10.0);
391        assert!((val - expected).abs() < 1e-10);
392    }
393
394    #[test]
395    fn test_lorentzian_symmetry() {
396        let left = lorentzian(90.0, 100.0, 10.0);
397        let right = lorentzian(110.0, 100.0, 10.0);
398        assert!((left - right).abs() < 1e-10);
399    }
400
401    #[test]
402    fn test_atomic_masses_known_elements() {
403        assert!((atomic_mass(1) - 1.00794).abs() < 0.001);
404        assert!((atomic_mass(6) - 12.011).abs() < 0.01);
405        assert!((atomic_mass(8) - 15.999).abs() < 0.01);
406    }
407
408    #[test]
409    fn test_h2_vibrational_analysis() {
410        let elements = [1u8, 1];
411        let positions = [[0.0, 0.0, 0.0], [0.74, 0.0, 0.0]];
412        let analysis =
413            compute_vibrational_analysis(&elements, &positions, HessianMethod::Xtb, Some(0.005))
414                .unwrap();
415
416        assert_eq!(analysis.n_atoms, 2);
417        assert_eq!(analysis.modes.len(), 6); // 3N = 6
418
419        // H₂ is linear: 5 zero + 1 stretch
420        // Should have at least 1 real mode (H-H stretch)
421        assert!(
422            analysis.n_real_modes >= 1,
423            "H₂ should have at least 1 real vibrational mode, got {}",
424            analysis.n_real_modes
425        );
426
427        // The stretch frequency should be positive
428        let real_modes: Vec<&VibrationalMode> =
429            analysis.modes.iter().filter(|m| m.is_real).collect();
430        assert!(!real_modes.is_empty(), "Should have at least one real mode");
431
432        // Check that ZPVE is positive for real modes
433        if analysis.n_real_modes > 0 {
434            assert!(analysis.zpve_ev > 0.0, "ZPVE should be positive");
435        }
436    }
437
438    #[test]
439    fn test_water_vibrational_analysis() {
440        let elements = [8u8, 1, 1];
441        let positions = [[0.0, 0.0, 0.0], [0.757, 0.586, 0.0], [-0.757, 0.586, 0.0]];
442        let analysis =
443            compute_vibrational_analysis(&elements, &positions, HessianMethod::Xtb, Some(0.005))
444                .unwrap();
445
446        assert_eq!(analysis.n_atoms, 3);
447        assert_eq!(analysis.modes.len(), 9); // 3N = 9
448
449        // Water is nonlinear: 6 zero + 3 real modes (sym stretch, asym stretch, bend)
450        // Due to numerical noise some translational modes may appear as real
451        // so we just check that we have at least 3 real modes
452        assert!(
453            analysis.n_real_modes >= 2,
454            "Water should have at least 2-3 real modes, got {}",
455            analysis.n_real_modes
456        );
457    }
458
459    #[test]
460    fn test_ir_spectrum_generation() {
461        let elements = [8u8, 1, 1];
462        let positions = [[0.0, 0.0, 0.0], [0.757, 0.586, 0.0], [-0.757, 0.586, 0.0]];
463        let analysis =
464            compute_vibrational_analysis(&elements, &positions, HessianMethod::Xtb, Some(0.005))
465                .unwrap();
466
467        let spectrum = compute_ir_spectrum(&analysis, 20.0, 400.0, 4000.0, 500);
468
469        assert_eq!(spectrum.wavenumbers.len(), 500);
470        assert_eq!(spectrum.intensities.len(), 500);
471        assert!(!spectrum.peaks.is_empty(), "Should have IR peaks");
472        assert!(
473            spectrum.intensities.iter().any(|&i| i > 0.0),
474            "Spectrum should have non-zero intensity"
475        );
476    }
477}