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    /// Thermochemistry at 298.15 K (if computed).
70    pub thermochemistry: Option<Thermochemistry>,
71    /// Semiempirical method used for Hessian.
72    pub method: String,
73    /// Notes and caveats.
74    pub notes: Vec<String>,
75}
76
77/// Thermochemical properties from RRHO approximation.
78#[derive(Debug, Clone, Serialize, Deserialize)]
79pub struct Thermochemistry {
80    /// Temperature in K.
81    pub temperature_k: f64,
82    /// Zero-point energy in kcal/mol.
83    pub zpve_kcal: f64,
84    /// Thermal energy correction in kcal/mol.
85    pub thermal_energy_kcal: f64,
86    /// Thermal enthalpy correction in kcal/mol (E_thermal + RT).
87    pub enthalpy_correction_kcal: f64,
88    /// Vibrational entropy in cal/(mol·K).
89    pub entropy_vib_cal: f64,
90    /// Gibbs free energy correction in kcal/mol (H - TS).
91    pub gibbs_correction_kcal: f64,
92}
93
94/// A peak in the IR spectrum.
95#[derive(Debug, Clone, Serialize, Deserialize)]
96pub struct IrPeak {
97    /// Frequency in cm⁻¹.
98    pub frequency_cm1: f64,
99    /// IR intensity in km/mol.
100    pub ir_intensity: f64,
101    /// Mode index.
102    pub mode_index: usize,
103    /// Functional group assignment (e.g., "O-H stretch", "C=O stretch").
104    pub assignment: String,
105}
106
107/// Broadened IR spectrum.
108#[derive(Debug, Clone, Serialize, Deserialize)]
109pub struct IrSpectrum {
110    /// Wavenumber axis (cm⁻¹).
111    pub wavenumbers: Vec<f64>,
112    /// Transmittance/absorbance axis.
113    pub intensities: Vec<f64>,
114    /// Discrete peaks from vibrational analysis.
115    pub peaks: Vec<IrPeak>,
116    /// Broadening width in cm⁻¹.
117    pub gamma: f64,
118    /// Notes.
119    pub notes: Vec<String>,
120}
121
122// Constants for unit conversion
123// 1 eV = 8065.544 cm⁻¹
124const EV_TO_CM1: f64 = 8065.544;
125// 1 amu = 1.66054e-27 kg, 1 eV = 1.60218e-19 J, 1 Å = 1e-10 m
126// Hessian in eV/Ų, mass in amu → frequency in cm⁻¹:
127// ω² = H/(m) in eV/(amu·Å²) → convert to s⁻²:
128// factor = eV_to_J / (amu_to_kg * Å_to_m²) = 1.60218e-19 / (1.66054e-27 * 1e-20) = 9.6485e27
129// ν_cm1 = sqrt(factor * eigenvalue) / (2π * c) where c = 2.998e10 cm/s
130const HESSIAN_TO_FREQUENCY_FACTOR: f64 = 9.6485e27; // eV/(amu·Å²) → s⁻²
131const INV_2PI_C: f64 = 1.0 / (2.0 * std::f64::consts::PI * 2.99792458e10); // 1/(2πc) in s/cm
132
133/// Compute the molecular dipole for use in IR intensity calculation.
134fn compute_dipole_vector(
135    elements: &[u8],
136    positions: &[[f64; 3]],
137    method: HessianMethod,
138) -> Result<[f64; 3], String> {
139    match method {
140        HessianMethod::Eht => {
141            let eht = crate::eht::solve_eht(elements, positions, None)?;
142            let dipole = crate::dipole::compute_dipole_from_eht(
143                elements,
144                positions,
145                &eht.coefficients,
146                eht.n_electrons,
147            );
148            Ok(dipole.vector)
149        }
150        HessianMethod::Pm3 => {
151            let pm3 = crate::pm3::solve_pm3(elements, positions)?;
152            // Use Mulliken charges for dipole
153            let dipole = crate::dipole::compute_dipole(&pm3.mulliken_charges, positions);
154            Ok(dipole.vector)
155        }
156        HessianMethod::Xtb => {
157            let xtb = crate::xtb::solve_xtb(elements, positions)?;
158            let dipole = crate::dipole::compute_dipole(&xtb.mulliken_charges, positions);
159            Ok(dipole.vector)
160        }
161        HessianMethod::Uff => {
162            // UFF has no electronic structure; use Gasteiger charges as approximation
163            let n = elements.len();
164            let charges = vec![0.0; n]; // neutral approximation
165            let dipole = crate::dipole::compute_dipole(&charges, positions);
166            Ok(dipole.vector)
167        }
168    }
169}
170
171/// Compute IR intensities from numerical dipole derivatives along normal modes.
172///
173/// I_k ∝ |dμ/dQ_k|² where Q_k is the k-th normal coordinate.
174fn compute_ir_intensities(
175    elements: &[u8],
176    positions: &[[f64; 3]],
177    normal_modes: &DMatrix<f64>,
178    method: HessianMethod,
179    delta: f64,
180) -> Result<Vec<f64>, String> {
181    let n_atoms = elements.len();
182    let n_modes = normal_modes.ncols();
183    let masses: Vec<f64> = elements.iter().map(|&z| atomic_mass(z)).collect();
184
185    let mut intensities = Vec::with_capacity(n_modes);
186
187    for k in 0..n_modes {
188        // Displace along normal mode k
189        let mut pos_plus: Vec<[f64; 3]> = positions.to_vec();
190        let mut pos_minus: Vec<[f64; 3]> = positions.to_vec();
191
192        for i in 0..n_atoms {
193            let sqrt_m = masses[i].sqrt();
194            for c in 0..3 {
195                let disp = normal_modes[(3 * i + c, k)] * delta / sqrt_m;
196                pos_plus[i][c] += disp;
197                pos_minus[i][c] -= disp;
198            }
199        }
200
201        let mu_plus = compute_dipole_vector(elements, &pos_plus, method)?;
202        let mu_minus = compute_dipole_vector(elements, &pos_minus, method)?;
203
204        // dμ/dQ ≈ (μ+ - μ-) / (2δ)
205        let dmu_dq: Vec<f64> = (0..3)
206            .map(|c| (mu_plus[c] - mu_minus[c]) / (2.0 * delta))
207            .collect();
208
209        // I ∝ |dμ/dQ|² — convert to km/mol with empirical scaling
210        let intensity = dmu_dq.iter().map(|x| x * x).sum::<f64>();
211        // Scale: Debye²/Ų·amu → km/mol (empirical: ~42.256)
212        intensities.push(intensity * 42.256);
213    }
214
215    Ok(intensities)
216}
217
218/// Lorentzian line shape: L(x) = (1/π) · γ / [(x - x₀)² + γ²]
219fn lorentzian(x: f64, x0: f64, gamma: f64) -> f64 {
220    let g = gamma.max(1e-6);
221    let dx = x - x0;
222    g / (std::f64::consts::PI * (dx * dx + g * g))
223}
224
225/// Gaussian line shape: G(x) = (1/(σ√(2π))) · exp(-(x-x₀)²/(2σ²))
226fn gaussian(x: f64, x0: f64, sigma: f64) -> f64 {
227    let s = sigma.max(1e-6);
228    let dx = x - x0;
229    (-(dx * dx) / (2.0 * s * s)).exp() / (s * (2.0 * std::f64::consts::PI).sqrt())
230}
231
232/// Assign functional group label to an IR frequency based on standard ranges.
233fn assign_ir_peak(frequency_cm1: f64) -> String {
234    if frequency_cm1 > 3500.0 {
235        "O-H stretch (broad)".to_string()
236    } else if frequency_cm1 > 3300.0 {
237        "N-H stretch".to_string()
238    } else if frequency_cm1 > 3000.0 {
239        "sp2 C-H stretch".to_string()
240    } else if frequency_cm1 > 2800.0 {
241        "sp3 C-H stretch".to_string()
242    } else if frequency_cm1 > 2100.0 && frequency_cm1 < 2300.0 {
243        "C≡N or C≡C stretch".to_string()
244    } else if frequency_cm1 > 1680.0 && frequency_cm1 < 1800.0 {
245        "C=O stretch".to_string()
246    } else if frequency_cm1 > 1600.0 && frequency_cm1 < 1680.0 {
247        "C=C stretch".to_string()
248    } else if frequency_cm1 > 1400.0 && frequency_cm1 < 1600.0 {
249        "aromatic C=C stretch".to_string()
250    } else if frequency_cm1 > 1300.0 && frequency_cm1 < 1400.0 {
251        "C-H bend".to_string()
252    } else if frequency_cm1 > 1000.0 && frequency_cm1 < 1300.0 {
253        "C-O stretch".to_string()
254    } else if frequency_cm1 > 600.0 && frequency_cm1 < 900.0 {
255        "C-H out-of-plane bend".to_string()
256    } else {
257        "skeletal mode".to_string()
258    }
259}
260
261/// Compute RRHO thermochemistry from vibrational frequencies.
262///
263/// Uses the rigid-rotor harmonic-oscillator approximation at the given temperature.
264fn compute_thermochemistry(frequencies_cm1: &[f64], temperature_k: f64) -> Thermochemistry {
265    // Constants
266    const CM1_TO_EV: f64 = 1.0 / 8065.544;
267    const EV_TO_KCAL: f64 = 23.0605;
268    const R_KCAL: f64 = 0.001987204; // kcal/(mol·K)
269    const R_CAL: f64 = 1.987204; // cal/(mol·K)
270    const KB_EV: f64 = 8.617333e-5; // eV/K
271
272    let kbt_ev = KB_EV * temperature_k;
273
274    let real_freqs: Vec<f64> = frequencies_cm1
275        .iter()
276        .filter(|&&f| f > 50.0)
277        .copied()
278        .collect();
279
280    // ZPVE = Σ (1/2)hν with 0.9 anharmonic scaling
281    let zpve_ev: f64 = real_freqs.iter().map(|&f| 0.5 * f * CM1_TO_EV * 0.9).sum();
282    let zpve_kcal = zpve_ev * EV_TO_KCAL;
283
284    // Thermal energy (vibrational contribution): Σ hν / (exp(hν/kT) - 1)
285    let thermal_e_ev: f64 = real_freqs
286        .iter()
287        .map(|&f| {
288            let hnu = f * CM1_TO_EV;
289            let x = hnu / kbt_ev;
290            if x > 100.0 {
291                0.0
292            } else {
293                hnu / (x.exp() - 1.0)
294            }
295        })
296        .sum();
297    let thermal_energy_kcal = thermal_e_ev * EV_TO_KCAL;
298
299    // Enthalpy correction: E_thermal + RT
300    let enthalpy_correction_kcal = thermal_energy_kcal + R_KCAL * temperature_k;
301
302    // Vibrational entropy: Σ [x/(exp(x)-1) - ln(1-exp(-x))]·R
303    let entropy_vib_cal: f64 = real_freqs
304        .iter()
305        .map(|&f| {
306            let hnu = f * CM1_TO_EV;
307            let x = hnu / kbt_ev;
308            if x > 100.0 {
309                0.0
310            } else {
311                let ex = x.exp();
312                R_CAL * (x / (ex - 1.0) - (1.0 - (-x).exp()).ln())
313            }
314        })
315        .sum();
316
317    // Gibbs correction: H - TS
318    let gibbs_correction_kcal = enthalpy_correction_kcal - temperature_k * entropy_vib_cal / 1000.0;
319
320    Thermochemistry {
321        temperature_k,
322        zpve_kcal,
323        thermal_energy_kcal,
324        enthalpy_correction_kcal,
325        entropy_vib_cal,
326        gibbs_correction_kcal,
327    }
328}
329
330/// Perform a complete vibrational analysis from elements and positions.
331///
332/// Computes the numerical Hessian, mass-weighted eigenvalue problem,
333/// vibrational frequencies, and IR intensities.
334pub fn compute_vibrational_analysis(
335    elements: &[u8],
336    positions: &[[f64; 3]],
337    method: HessianMethod,
338    step_size: Option<f64>,
339) -> Result<VibrationalAnalysis, String> {
340    if method == HessianMethod::Uff {
341        return Err(
342            "UFF requires SMILES; use compute_vibrational_analysis_uff instead".to_string(),
343        );
344    }
345
346    let n_atoms = elements.len();
347    let delta = step_size.unwrap_or(0.005);
348
349    if n_atoms < 2 {
350        return Err("Need at least 2 atoms for vibrational analysis".to_string());
351    }
352
353    let hessian = compute_numerical_hessian(elements, positions, method, Some(delta))?;
354    build_vibrational_analysis_from_hessian(elements, positions, &hessian, method, delta)
355}
356
357/// Perform vibrational analysis using UFF analytical Hessian (gradient-difference method).
358///
359/// This avoids the expensive O(9N²) energy evaluations of the standard numerical Hessian
360/// by using O(6N) gradient evaluations with UFF's analytical gradients.
361pub fn compute_vibrational_analysis_uff(
362    smiles: &str,
363    elements: &[u8],
364    positions: &[[f64; 3]],
365    step_size: Option<f64>,
366) -> Result<VibrationalAnalysis, String> {
367    let n_atoms = elements.len();
368    let delta = step_size.unwrap_or(0.005);
369
370    if n_atoms < 2 {
371        return Err("Need at least 2 atoms for vibrational analysis".to_string());
372    }
373
374    let coords_flat: Vec<f64> = positions.iter().flat_map(|p| p.iter().copied()).collect();
375    let hessian =
376        super::hessian::compute_uff_analytical_hessian(smiles, &coords_flat, Some(delta))?;
377    build_vibrational_analysis_from_hessian(
378        elements,
379        positions,
380        &hessian,
381        HessianMethod::Uff,
382        delta,
383    )
384}
385
386/// Build vibrational analysis from a pre-computed Hessian matrix.
387fn build_vibrational_analysis_from_hessian(
388    elements: &[u8],
389    positions: &[[f64; 3]],
390    hessian: &DMatrix<f64>,
391    method: HessianMethod,
392    delta: f64,
393) -> Result<VibrationalAnalysis, String> {
394    let n_atoms = elements.len();
395    let n3 = 3 * n_atoms;
396
397    // 2. Build mass-weighted Hessian: H'_{ij} = H_{ij} / sqrt(m_i * m_j)
398    let masses: Vec<f64> = elements.iter().map(|&z| atomic_mass(z)).collect();
399    let mut mw_hessian = DMatrix::zeros(n3, n3);
400    for i in 0..n3 {
401        let mi = masses[i / 3];
402        for j in 0..n3 {
403            let mj = masses[j / 3];
404            mw_hessian[(i, j)] = hessian[(i, j)] / (mi * mj).sqrt();
405        }
406    }
407
408    // 3. Diagonalize mass-weighted Hessian
409    let eigen = mw_hessian.symmetric_eigen();
410
411    // 4. Sort eigenvalues and eigenvectors by eigenvalue
412    let mut indices: Vec<usize> = (0..n3).collect();
413    indices.sort_by(|&a, &b| {
414        eigen.eigenvalues[a]
415            .partial_cmp(&eigen.eigenvalues[b])
416            .unwrap_or(std::cmp::Ordering::Equal)
417    });
418
419    // 5. Convert eigenvalues to frequencies (cm⁻¹)
420    // For linear molecules: 5 zero modes, for nonlinear: 6
421    let is_linear = n_atoms == 2; // simplified check
422    let _n_tr = if is_linear { 5 } else { 6 };
423    let freq_threshold = 50.0; // cm⁻¹ threshold for "real" modes
424
425    let mut sorted_eigenvalues = Vec::with_capacity(n3);
426    let mut sorted_modes = DMatrix::zeros(n3, n3);
427    for (new_idx, &old_idx) in indices.iter().enumerate() {
428        sorted_eigenvalues.push(eigen.eigenvalues[old_idx]);
429        for i in 0..n3 {
430            sorted_modes[(i, new_idx)] = eigen.eigenvectors[(i, old_idx)];
431        }
432    }
433
434    // Convert eigenvalues to frequencies
435    let frequencies: Vec<f64> = sorted_eigenvalues
436        .iter()
437        .map(|&ev| {
438            if ev >= 0.0 {
439                // ν = (1/2πc) * sqrt(eigenvalue * factor)
440                (ev * HESSIAN_TO_FREQUENCY_FACTOR).sqrt() * INV_2PI_C
441            } else {
442                // Imaginary frequency (negative eigenvalue)
443                -((-ev) * HESSIAN_TO_FREQUENCY_FACTOR).sqrt() * INV_2PI_C
444            }
445        })
446        .collect();
447
448    // 6. Compute IR intensities from dipole derivatives
449    let ir_intensities = compute_ir_intensities(
450        elements,
451        positions,
452        &sorted_modes,
453        method,
454        delta * 2.0, // slightly larger step for dipole derivatives
455    )?;
456
457    // 7. Build mode list
458    let mut modes = Vec::with_capacity(n3);
459    let mut n_real = 0;
460    let mut zpve = 0.0;
461
462    for k in 0..n3 {
463        let freq = frequencies[k];
464        let is_real = freq.abs() > freq_threshold;
465        if is_real && freq > 0.0 {
466            n_real += 1;
467            // ZPVE contribution: (1/2)hν in eV
468            zpve += 0.5 * freq / EV_TO_CM1;
469        }
470
471        let displacement: Vec<f64> = (0..n3).map(|i| sorted_modes[(i, k)]).collect();
472
473        modes.push(VibrationalMode {
474            frequency_cm1: freq,
475            ir_intensity: ir_intensities.get(k).copied().unwrap_or(0.0),
476            displacement,
477            is_real,
478        });
479    }
480
481    let method_name = match method {
482        HessianMethod::Eht => "EHT",
483        HessianMethod::Pm3 => "PM3",
484        HessianMethod::Xtb => "xTB",
485        HessianMethod::Uff => "UFF",
486    };
487
488    // Compute thermochemistry at 298.15 K
489    let thermochemistry = Some(compute_thermochemistry(&frequencies, 298.15));
490
491    Ok(VibrationalAnalysis {
492        n_atoms,
493        modes,
494        n_real_modes: n_real,
495        zpve_ev: zpve,
496        thermochemistry,
497        method: method_name.to_string(),
498        notes: vec![
499            format!(
500                "Numerical Hessian computed with {} using central finite differences (δ = {} Å).",
501                method_name, delta
502            ),
503            "IR intensities derived from numerical dipole derivatives along normal coordinates."
504                .to_string(),
505            "Frequencies below 50 cm⁻¹ are classified as translation/rotation modes.".to_string(),
506        ],
507    })
508}
509
510/// Broadening function type for IR spectrum generation.
511#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
512pub enum BroadeningType {
513    /// Lorentzian line shape (default).
514    Lorentzian,
515    /// Gaussian line shape.
516    Gaussian,
517}
518
519/// Generate a broadened IR spectrum from vibrational analysis.
520///
521/// `gamma`: broadening width in cm⁻¹ (typically 10–30).
522/// `wn_min`, `wn_max`: wavenumber range in cm⁻¹.
523/// `n_points`: number of grid points.
524pub fn compute_ir_spectrum(
525    analysis: &VibrationalAnalysis,
526    gamma: f64,
527    wn_min: f64,
528    wn_max: f64,
529    n_points: usize,
530) -> IrSpectrum {
531    compute_ir_spectrum_with_broadening(
532        analysis,
533        gamma,
534        wn_min,
535        wn_max,
536        n_points,
537        BroadeningType::Lorentzian,
538    )
539}
540
541/// Generate a broadened IR spectrum with selectable broadening function.
542pub fn compute_ir_spectrum_with_broadening(
543    analysis: &VibrationalAnalysis,
544    gamma: f64,
545    wn_min: f64,
546    wn_max: f64,
547    n_points: usize,
548    broadening: BroadeningType,
549) -> IrSpectrum {
550    let n_points = n_points.max(2);
551    let step = (wn_max - wn_min) / (n_points as f64 - 1.0);
552    let wavenumbers: Vec<f64> = (0..n_points).map(|i| wn_min + step * i as f64).collect();
553    let mut intensities = vec![0.0; n_points];
554
555    let mut peaks = Vec::new();
556
557    for (mode_idx, mode) in analysis.modes.iter().enumerate() {
558        if !mode.is_real || mode.frequency_cm1 <= 0.0 {
559            continue;
560        }
561
562        peaks.push(IrPeak {
563            frequency_cm1: mode.frequency_cm1,
564            ir_intensity: mode.ir_intensity,
565            mode_index: mode_idx,
566            assignment: assign_ir_peak(mode.frequency_cm1),
567        });
568
569        for (idx, &wn) in wavenumbers.iter().enumerate() {
570            intensities[idx] += mode.ir_intensity
571                * match broadening {
572                    BroadeningType::Lorentzian => lorentzian(wn, mode.frequency_cm1, gamma),
573                    BroadeningType::Gaussian => gaussian(wn, mode.frequency_cm1, gamma),
574                };
575        }
576    }
577
578    peaks.sort_by(|a, b| {
579        b.ir_intensity
580            .partial_cmp(&a.ir_intensity)
581            .unwrap_or(std::cmp::Ordering::Equal)
582    });
583
584    IrSpectrum {
585        wavenumbers,
586        intensities,
587        peaks,
588        gamma,
589        notes: vec![
590            format!(
591                "IR spectrum generated with Lorentzian broadening (γ = {} cm⁻¹).",
592                gamma
593            ),
594            format!("Vibrational analysis method: {}.", analysis.method),
595        ],
596    }
597}
598
599#[cfg(test)]
600mod tests {
601    use super::*;
602
603    #[test]
604    fn test_lorentzian_peak_at_center() {
605        let val = lorentzian(100.0, 100.0, 10.0);
606        // At center: L = 1/(π·γ)
607        let expected = 1.0 / (std::f64::consts::PI * 10.0);
608        assert!((val - expected).abs() < 1e-10);
609    }
610
611    #[test]
612    fn test_lorentzian_symmetry() {
613        let left = lorentzian(90.0, 100.0, 10.0);
614        let right = lorentzian(110.0, 100.0, 10.0);
615        assert!((left - right).abs() < 1e-10);
616    }
617
618    #[test]
619    fn test_atomic_masses_known_elements() {
620        assert!((atomic_mass(1) - 1.00794).abs() < 0.001);
621        assert!((atomic_mass(6) - 12.011).abs() < 0.01);
622        assert!((atomic_mass(8) - 15.999).abs() < 0.01);
623    }
624
625    #[test]
626    fn test_h2_vibrational_analysis() {
627        let elements = [1u8, 1];
628        let positions = [[0.0, 0.0, 0.0], [0.74, 0.0, 0.0]];
629        let analysis =
630            compute_vibrational_analysis(&elements, &positions, HessianMethod::Xtb, Some(0.005))
631                .unwrap();
632
633        assert_eq!(analysis.n_atoms, 2);
634        assert_eq!(analysis.modes.len(), 6); // 3N = 6
635
636        // H₂ is linear: 5 zero + 1 stretch
637        // Should have at least 1 real mode (H-H stretch)
638        assert!(
639            analysis.n_real_modes >= 1,
640            "H₂ should have at least 1 real vibrational mode, got {}",
641            analysis.n_real_modes
642        );
643
644        // The stretch frequency should be positive
645        let real_modes: Vec<&VibrationalMode> =
646            analysis.modes.iter().filter(|m| m.is_real).collect();
647        assert!(!real_modes.is_empty(), "Should have at least one real mode");
648
649        // Check that ZPVE is positive for real modes
650        if analysis.n_real_modes > 0 {
651            assert!(analysis.zpve_ev > 0.0, "ZPVE should be positive");
652        }
653    }
654
655    #[test]
656    fn test_water_vibrational_analysis() {
657        let elements = [8u8, 1, 1];
658        let positions = [[0.0, 0.0, 0.0], [0.757, 0.586, 0.0], [-0.757, 0.586, 0.0]];
659        let analysis =
660            compute_vibrational_analysis(&elements, &positions, HessianMethod::Xtb, Some(0.005))
661                .unwrap();
662
663        assert_eq!(analysis.n_atoms, 3);
664        assert_eq!(analysis.modes.len(), 9); // 3N = 9
665
666        // Water is nonlinear: 6 zero + 3 real modes (sym stretch, asym stretch, bend)
667        // Due to numerical noise some translational modes may appear as real
668        // so we just check that we have at least 3 real modes
669        assert!(
670            analysis.n_real_modes >= 2,
671            "Water should have at least 2-3 real modes, got {}",
672            analysis.n_real_modes
673        );
674    }
675
676    #[test]
677    fn test_ir_spectrum_generation() {
678        let elements = [8u8, 1, 1];
679        let positions = [[0.0, 0.0, 0.0], [0.757, 0.586, 0.0], [-0.757, 0.586, 0.0]];
680        let analysis =
681            compute_vibrational_analysis(&elements, &positions, HessianMethod::Xtb, Some(0.005))
682                .unwrap();
683
684        let spectrum = compute_ir_spectrum(&analysis, 20.0, 400.0, 4000.0, 500);
685
686        assert_eq!(spectrum.wavenumbers.len(), 500);
687        assert_eq!(spectrum.intensities.len(), 500);
688        assert!(!spectrum.peaks.is_empty(), "Should have IR peaks");
689        assert!(
690            spectrum.intensities.iter().any(|&i| i > 0.0),
691            "Spectrum should have non-zero intensity"
692        );
693    }
694}