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        3 => 6.941,
18        4 => 9.01218,
19        5 => 10.811,
20        6 => 12.0107,
21        7 => 14.0067,
22        8 => 15.9994,
23        9 => 18.9984,
24        10 => 20.1797,
25        11 => 22.9898,
26        12 => 24.305,
27        13 => 26.9815,
28        14 => 28.0855,
29        15 => 30.9738,
30        16 => 32.065,
31        17 => 35.453,
32        18 => 39.948,
33        19 => 39.0983,
34        20 => 40.078,
35        21 => 44.9559,
36        22 => 47.867,
37        23 => 50.9415,
38        24 => 51.9961,
39        25 => 54.9380,
40        26 => 55.845,
41        27 => 58.9332,
42        28 => 58.6934,
43        29 => 63.546,
44        30 => 65.38,
45        31 => 69.723,
46        32 => 72.630,
47        33 => 74.9216,
48        34 => 78.971,
49        35 => 79.904,
50        36 => 83.798,
51        37 => 85.4678,
52        38 => 87.62,
53        39 => 88.9058,
54        40 => 91.224,
55        41 => 92.9064,
56        42 => 95.95,
57        43 => 98.0,
58        44 => 101.07,
59        45 => 102.906,
60        46 => 106.42,
61        47 => 107.868,
62        48 => 112.414,
63        49 => 114.818,
64        50 => 118.710,
65        51 => 121.760,
66        52 => 127.60,
67        53 => 126.904,
68        54 => 131.293,
69        55 => 132.905,
70        56 => 137.327,
71        72 => 178.49,
72        73 => 180.948,
73        74 => 183.84,
74        75 => 186.207,
75        76 => 190.23,
76        77 => 192.217,
77        78 => 195.084,
78        79 => 196.967,
79        80 => 200.592,
80        81 => 204.38,
81        82 => 207.2,
82        83 => 208.980,
83        // Lanthanides (La-Lu)
84        57 => 138.905,
85        58 => 140.116,
86        59 => 140.908,
87        60 => 144.242,
88        61 => 145.0,
89        62 => 150.36,
90        63 => 151.964,
91        64 => 157.25,
92        65 => 158.925,
93        66 => 162.500,
94        67 => 164.930,
95        68 => 167.259,
96        69 => 168.934,
97        70 => 173.045,
98        71 => 174.967,
99        // Remaining heavy elements
100        84 => 209.0,
101        85 => 210.0,
102        86 => 222.0,
103        _ => z as f64 * 1.5, // rough fallback for actinides and transactinides
104    }
105}
106
107/// A single vibrational normal mode.
108#[derive(Debug, Clone, Serialize, Deserialize)]
109pub struct VibrationalMode {
110    /// Frequency in cm⁻¹ (negative = imaginary).
111    pub frequency_cm1: f64,
112    /// IR intensity in km/mol.
113    pub ir_intensity: f64,
114    /// Displacement vector (3N elements, mass-weighted normal coordinate).
115    pub displacement: Vec<f64>,
116    /// Whether this is a real mode (true) or translation/rotation (false).
117    pub is_real: bool,
118}
119
120/// Complete vibrational analysis result.
121#[derive(Debug, Clone, Serialize, Deserialize)]
122pub struct VibrationalAnalysis {
123    /// Number of atoms.
124    pub n_atoms: usize,
125    /// Atomic numbers carried through to improve downstream peak assignment.
126    #[serde(default)]
127    pub elements: Vec<u8>,
128    /// All vibrational modes sorted by frequency.
129    pub modes: Vec<VibrationalMode>,
130    /// Number of real vibrational modes (excluding translation/rotation).
131    pub n_real_modes: usize,
132    /// Zero-point vibrational energy in eV.
133    pub zpve_ev: f64,
134    /// Thermochemistry at 298.15 K (if computed).
135    pub thermochemistry: Option<Thermochemistry>,
136    /// Semiempirical method used for Hessian.
137    pub method: String,
138    /// Notes and caveats.
139    pub notes: Vec<String>,
140}
141
142/// Thermochemical properties from RRHO approximation.
143#[derive(Debug, Clone, Serialize, Deserialize)]
144pub struct Thermochemistry {
145    /// Temperature in K.
146    pub temperature_k: f64,
147    /// Zero-point energy in kcal/mol.
148    pub zpve_kcal: f64,
149    /// Thermal energy correction in kcal/mol.
150    pub thermal_energy_kcal: f64,
151    /// Thermal enthalpy correction in kcal/mol (E_thermal + RT).
152    pub enthalpy_correction_kcal: f64,
153    /// Vibrational entropy in cal/(mol·K).
154    pub entropy_vib_cal: f64,
155    /// Gibbs free energy correction in kcal/mol (H - TS).
156    pub gibbs_correction_kcal: f64,
157}
158
159/// A peak in the IR spectrum.
160#[derive(Debug, Clone, Serialize, Deserialize)]
161pub struct IrPeak {
162    /// Frequency in cm⁻¹.
163    pub frequency_cm1: f64,
164    /// IR intensity in km/mol.
165    pub ir_intensity: f64,
166    /// Mode index.
167    pub mode_index: usize,
168    /// Functional group assignment (e.g., "O-H stretch", "C=O stretch").
169    pub assignment: String,
170}
171
172/// Broadened IR spectrum.
173#[derive(Debug, Clone, Serialize, Deserialize)]
174pub struct IrSpectrum {
175    /// Wavenumber axis (cm⁻¹).
176    pub wavenumbers: Vec<f64>,
177    /// Transmittance/absorbance axis.
178    pub intensities: Vec<f64>,
179    /// Discrete peaks from vibrational analysis.
180    pub peaks: Vec<IrPeak>,
181    /// Broadening width in cm⁻¹.
182    pub gamma: f64,
183    /// Notes.
184    pub notes: Vec<String>,
185}
186
187// Constants for unit conversion
188// 1 eV = 8065.544 cm⁻¹
189const EV_TO_CM1: f64 = 8065.544;
190const KCAL_MOL_PER_EV: f64 = 23.0605;
191// 1 amu = 1.66054e-27 kg, 1 eV = 1.60218e-19 J, 1 Å = 1e-10 m
192// Hessian in eV/Ų, mass in amu → frequency in cm⁻¹:
193// ω² = H/(m) in eV/(amu·Å²) → convert to s⁻²:
194// factor = eV_to_J / (amu_to_kg * Å_to_m²) = 1.60218e-19 / (1.66054e-27 * 1e-20) = 9.6485e27
195// ν_cm1 = sqrt(factor * eigenvalue) / (2π * c) where c = 2.998e10 cm/s
196const HESSIAN_TO_FREQUENCY_FACTOR: f64 = 9.6485e27; // eV/(amu·Å²) → s⁻²
197const INV_2PI_C: f64 = 1.0 / (2.0 * std::f64::consts::PI * 2.99792458e10); // 1/(2πc) in s/cm
198
199fn hessian_frequency_factor(method: HessianMethod) -> f64 {
200    match method {
201        HessianMethod::Uff => HESSIAN_TO_FREQUENCY_FACTOR / KCAL_MOL_PER_EV,
202        _ => HESSIAN_TO_FREQUENCY_FACTOR,
203    }
204}
205
206fn push_orthonormal_basis_vector(basis: &mut Vec<Vec<f64>>, mut candidate: Vec<f64>) {
207    for existing in basis.iter() {
208        let dot = candidate
209            .iter()
210            .zip(existing.iter())
211            .map(|(left, right)| left * right)
212            .sum::<f64>();
213
214        for (value, existing_value) in candidate.iter_mut().zip(existing.iter()) {
215            *value -= dot * existing_value;
216        }
217    }
218
219    let norm = candidate
220        .iter()
221        .map(|value| value * value)
222        .sum::<f64>()
223        .sqrt();
224    if norm <= 1e-8 {
225        return;
226    }
227
228    for value in &mut candidate {
229        *value /= norm;
230    }
231
232    basis.push(candidate);
233}
234
235fn build_rigid_body_basis(masses: &[f64], positions: &[[f64; 3]]) -> DMatrix<f64> {
236    let n_atoms = positions.len();
237    let n3 = 3 * n_atoms;
238
239    let total_mass = masses.iter().sum::<f64>().max(1e-12);
240    let mut center_of_mass = [0.0; 3];
241    for (mass, position) in masses.iter().zip(positions.iter()) {
242        for axis in 0..3 {
243            center_of_mass[axis] += mass * position[axis];
244        }
245    }
246    for axis in 0..3 {
247        center_of_mass[axis] /= total_mass;
248    }
249
250    let centered_positions: Vec<[f64; 3]> = positions
251        .iter()
252        .map(|position| {
253            [
254                position[0] - center_of_mass[0],
255                position[1] - center_of_mass[1],
256                position[2] - center_of_mass[2],
257            ]
258        })
259        .collect();
260
261    let mut basis_vectors = Vec::with_capacity(6);
262
263    for axis in 0..3 {
264        let mut translation = vec![0.0; n3];
265        for (atom_index, mass) in masses.iter().enumerate() {
266            translation[3 * atom_index + axis] = mass.sqrt();
267        }
268        push_orthonormal_basis_vector(&mut basis_vectors, translation);
269    }
270
271    for axis in 0..3 {
272        let mut rotation = vec![0.0; n3];
273        for (atom_index, (mass, position)) in
274            masses.iter().zip(centered_positions.iter()).enumerate()
275        {
276            let scaled = mass.sqrt();
277            let components = match axis {
278                0 => [0.0, -position[2], position[1]],
279                1 => [position[2], 0.0, -position[0]],
280                _ => [-position[1], position[0], 0.0],
281            };
282
283            for component in 0..3 {
284                rotation[3 * atom_index + component] = scaled * components[component];
285            }
286        }
287        push_orthonormal_basis_vector(&mut basis_vectors, rotation);
288    }
289
290    let mut basis = DMatrix::zeros(n3, basis_vectors.len());
291    for (column, vector) in basis_vectors.iter().enumerate() {
292        for (row, value) in vector.iter().enumerate() {
293            basis[(row, column)] = *value;
294        }
295    }
296
297    basis
298}
299
300fn project_rigid_body_modes(mw_hessian: &DMatrix<f64>, rigid_basis: &DMatrix<f64>) -> DMatrix<f64> {
301    if rigid_basis.ncols() == 0 {
302        return mw_hessian.clone();
303    }
304
305    let n = mw_hessian.nrows();
306    let projector = DMatrix::<f64>::identity(n, n) - rigid_basis * rigid_basis.transpose();
307    let mut projected = &projector * mw_hessian * &projector;
308
309    for i in 0..n {
310        for j in (i + 1)..n {
311            let average = 0.5 * (projected[(i, j)] + projected[(j, i)]);
312            projected[(i, j)] = average;
313            projected[(j, i)] = average;
314        }
315    }
316
317    projected
318}
319
320fn rigid_body_overlap(displacement: &[f64], rigid_basis: &DMatrix<f64>) -> f64 {
321    if rigid_basis.ncols() == 0 {
322        return 0.0;
323    }
324
325    let mut overlap = 0.0;
326    for column in 0..rigid_basis.ncols() {
327        let mut dot = 0.0;
328        for row in 0..displacement.len() {
329            dot += displacement[row] * rigid_basis[(row, column)];
330        }
331        overlap += dot * dot;
332    }
333
334    overlap.clamp(0.0, 1.0)
335}
336
337fn relax_uff_geometry(smiles: &str, positions: &[[f64; 3]]) -> Result<Vec<[f64; 3]>, String> {
338    let mol = crate::graph::Molecule::from_smiles(smiles)?;
339    if mol.graph.node_count() != positions.len() {
340        return Err("SMILES atom count does not match coordinates for UFF relaxation".to_string());
341    }
342
343    let ff = crate::forcefield::builder::build_uff_force_field(&mol);
344    let mut coords: Vec<f64> = positions
345        .iter()
346        .flat_map(|position| position.iter().copied())
347        .collect();
348    let mut gradient = vec![0.0; coords.len()];
349    let mut energy = ff.compute_system_energy_and_gradients(&coords, &mut gradient);
350
351    for _ in 0..64 {
352        let grad_norm = gradient
353            .iter()
354            .map(|value| value * value)
355            .sum::<f64>()
356            .sqrt();
357        if grad_norm < 1e-4 {
358            break;
359        }
360
361        let max_component = gradient
362            .iter()
363            .map(|value| value.abs())
364            .fold(0.0_f64, f64::max);
365        let gradient_scale = if max_component > 10.0 {
366            10.0 / max_component
367        } else {
368            1.0
369        };
370
371        let mut step = 0.02;
372        let mut improved = false;
373
374        while step >= 1e-6 {
375            let trial_coords: Vec<f64> = coords
376                .iter()
377                .zip(gradient.iter())
378                .map(|(coord, grad)| coord - step * gradient_scale * grad)
379                .collect();
380            let mut trial_gradient = vec![0.0; coords.len()];
381            let trial_energy =
382                ff.compute_system_energy_and_gradients(&trial_coords, &mut trial_gradient);
383
384            if trial_energy < energy {
385                coords = trial_coords;
386                gradient = trial_gradient;
387                energy = trial_energy;
388                improved = true;
389                break;
390            }
391
392            step *= 0.5;
393        }
394
395        if !improved {
396            break;
397        }
398    }
399
400    Ok(coords
401        .chunks(3)
402        .map(|chunk| [chunk[0], chunk[1], chunk[2]])
403        .collect())
404}
405
406fn compute_uff_numerical_hessian(
407    smiles: &str,
408    coords_flat: &[f64],
409    delta: f64,
410) -> Result<DMatrix<f64>, String> {
411    let n3 = coords_flat.len();
412    let e0 = crate::compute_uff_energy(smiles, coords_flat)?;
413    let delta_sq = delta * delta;
414    let mut hessian = DMatrix::zeros(n3, n3);
415
416    for i in 0..n3 {
417        let mut coords_plus = coords_flat.to_vec();
418        let mut coords_minus = coords_flat.to_vec();
419        coords_plus[i] += delta;
420        coords_minus[i] -= delta;
421
422        let e_plus = crate::compute_uff_energy(smiles, &coords_plus)?;
423        let e_minus = crate::compute_uff_energy(smiles, &coords_minus)?;
424        hessian[(i, i)] = (e_plus - 2.0 * e0 + e_minus) / delta_sq;
425    }
426
427    for i in 0..n3 {
428        for j in (i + 1)..n3 {
429            let mut coords_pp = coords_flat.to_vec();
430            let mut coords_pm = coords_flat.to_vec();
431            let mut coords_mp = coords_flat.to_vec();
432            let mut coords_mm = coords_flat.to_vec();
433
434            coords_pp[i] += delta;
435            coords_pp[j] += delta;
436            coords_pm[i] += delta;
437            coords_pm[j] -= delta;
438            coords_mp[i] -= delta;
439            coords_mp[j] += delta;
440            coords_mm[i] -= delta;
441            coords_mm[j] -= delta;
442
443            let e_pp = crate::compute_uff_energy(smiles, &coords_pp)?;
444            let e_pm = crate::compute_uff_energy(smiles, &coords_pm)?;
445            let e_mp = crate::compute_uff_energy(smiles, &coords_mp)?;
446            let e_mm = crate::compute_uff_energy(smiles, &coords_mm)?;
447
448            let value = (e_pp - e_pm - e_mp + e_mm) / (4.0 * delta_sq);
449            hessian[(i, j)] = value;
450            hessian[(j, i)] = value;
451        }
452    }
453
454    for i in 0..n3 {
455        for j in (i + 1)..n3 {
456            let average = 0.5 * (hessian[(i, j)] + hessian[(j, i)]);
457            hessian[(i, j)] = average;
458            hessian[(j, i)] = average;
459        }
460    }
461
462    Ok(hessian)
463}
464
465fn has_significant_imaginary_vibrational_mode(analysis: &VibrationalAnalysis) -> bool {
466    analysis
467        .modes
468        .iter()
469        .any(|mode| mode.is_real && mode.frequency_cm1 < -100.0)
470}
471
472/// Compute the molecular dipole for use in IR intensity calculation.
473fn compute_dipole_vector(
474    elements: &[u8],
475    positions: &[[f64; 3]],
476    method: HessianMethod,
477    uff_charges: Option<&[f64]>,
478) -> Result<[f64; 3], String> {
479    match method {
480        HessianMethod::Eht => {
481            let eht = crate::eht::solve_eht(elements, positions, None)?;
482            let dipole = crate::dipole::compute_dipole_from_eht(
483                elements,
484                positions,
485                &eht.coefficients,
486                eht.n_electrons,
487            );
488            Ok(dipole.vector)
489        }
490        HessianMethod::Pm3 => {
491            let pm3 = crate::pm3::solve_pm3(elements, positions)?;
492            // Use Mulliken charges for dipole
493            let dipole = crate::dipole::compute_dipole(&pm3.mulliken_charges, positions);
494            Ok(dipole.vector)
495        }
496        HessianMethod::Xtb => {
497            let xtb = crate::xtb::solve_xtb(elements, positions)?;
498            let dipole = crate::dipole::compute_dipole(&xtb.mulliken_charges, positions);
499            Ok(dipole.vector)
500        }
501        HessianMethod::Uff => {
502            // UFF has no electronic structure; use Gasteiger charges as approximation
503            let charges = uff_charges
504                .map(|values| values.to_vec())
505                .unwrap_or_else(|| vec![0.0; elements.len()]);
506            let dipole = crate::dipole::compute_dipole(&charges, positions);
507            Ok(dipole.vector)
508        }
509    }
510}
511
512/// Compute IR intensities from numerical dipole derivatives along normal modes.
513///
514/// I_k ∝ |dμ/dQ_k|² where Q_k is the k-th normal coordinate.
515fn compute_ir_intensities(
516    elements: &[u8],
517    positions: &[[f64; 3]],
518    normal_modes: &DMatrix<f64>,
519    method: HessianMethod,
520    delta: f64,
521    uff_charges: Option<&[f64]>,
522) -> Result<Vec<f64>, String> {
523    #[cfg(feature = "parallel")]
524    {
525        compute_ir_intensities_parallel(
526            elements,
527            positions,
528            normal_modes,
529            method,
530            delta,
531            uff_charges,
532        )
533    }
534    #[cfg(not(feature = "parallel"))]
535    {
536        compute_ir_intensities_sequential(
537            elements,
538            positions,
539            normal_modes,
540            method,
541            delta,
542            uff_charges,
543        )
544    }
545}
546
547#[cfg(not(feature = "parallel"))]
548fn compute_ir_intensities_sequential(
549    elements: &[u8],
550    positions: &[[f64; 3]],
551    normal_modes: &DMatrix<f64>,
552    method: HessianMethod,
553    delta: f64,
554    uff_charges: Option<&[f64]>,
555) -> Result<Vec<f64>, String> {
556    let n_atoms = elements.len();
557    let n_modes = normal_modes.ncols();
558    let masses: Vec<f64> = elements.iter().map(|&z| atomic_mass(z)).collect();
559
560    let mut intensities = Vec::with_capacity(n_modes);
561
562    for k in 0..n_modes {
563        let mut pos_plus: Vec<[f64; 3]> = positions.to_vec();
564        let mut pos_minus: Vec<[f64; 3]> = positions.to_vec();
565
566        for i in 0..n_atoms {
567            let sqrt_m = masses[i].sqrt();
568            for c in 0..3 {
569                let disp = normal_modes[(3 * i + c, k)] * delta / sqrt_m;
570                pos_plus[i][c] += disp;
571                pos_minus[i][c] -= disp;
572            }
573        }
574
575        let mu_plus = compute_dipole_vector(elements, &pos_plus, method, uff_charges)?;
576        let mu_minus = compute_dipole_vector(elements, &pos_minus, method, uff_charges)?;
577
578        let dmu_dq: Vec<f64> = (0..3)
579            .map(|c| (mu_plus[c] - mu_minus[c]) / (2.0 * delta))
580            .collect();
581
582        let intensity = dmu_dq.iter().map(|x| x * x).sum::<f64>();
583        intensities.push(intensity * 42.256);
584    }
585
586    Ok(intensities)
587}
588
589/// Parallel IR intensities via rayon — each mode's dipole derivative
590/// displacement is independent.
591#[cfg(feature = "parallel")]
592fn compute_ir_intensities_parallel(
593    elements: &[u8],
594    positions: &[[f64; 3]],
595    normal_modes: &DMatrix<f64>,
596    method: HessianMethod,
597    delta: f64,
598    uff_charges: Option<&[f64]>,
599) -> Result<Vec<f64>, String> {
600    use rayon::prelude::*;
601
602    let n_atoms = elements.len();
603    let n_modes = normal_modes.ncols();
604    let masses: Vec<f64> = elements.iter().map(|&z| atomic_mass(z)).collect();
605
606    let results: Vec<Result<f64, String>> = (0..n_modes)
607        .into_par_iter()
608        .map(|k| {
609            let mut pos_plus: Vec<[f64; 3]> = positions.to_vec();
610            let mut pos_minus: Vec<[f64; 3]> = positions.to_vec();
611
612            for i in 0..n_atoms {
613                let sqrt_m = masses[i].sqrt();
614                for c in 0..3 {
615                    let disp = normal_modes[(3 * i + c, k)] * delta / sqrt_m;
616                    pos_plus[i][c] += disp;
617                    pos_minus[i][c] -= disp;
618                }
619            }
620
621            let mu_plus = compute_dipole_vector(elements, &pos_plus, method, uff_charges)?;
622            let mu_minus = compute_dipole_vector(elements, &pos_minus, method, uff_charges)?;
623
624            let dmu_dq: Vec<f64> = (0..3)
625                .map(|c| (mu_plus[c] - mu_minus[c]) / (2.0 * delta))
626                .collect();
627
628            let intensity = dmu_dq.iter().map(|x| x * x).sum::<f64>();
629            Ok(intensity * 42.256)
630        })
631        .collect();
632
633    results.into_iter().collect()
634}
635
636/// Lorentzian line shape: L(x) = (1/π) · γ / [(x - x₀)² + γ²]
637fn lorentzian(x: f64, x0: f64, gamma: f64) -> f64 {
638    let g = gamma.max(1e-6);
639    let dx = x - x0;
640    g / (std::f64::consts::PI * (dx * dx + g * g))
641}
642
643/// Gaussian line shape: G(x) = (1/(σ√(2π))) · exp(-(x-x₀)²/(2σ²))
644fn gaussian(x: f64, x0: f64, sigma: f64) -> f64 {
645    let s = sigma.max(1e-6);
646    let dx = x - x0;
647    (-(dx * dx) / (2.0 * s * s)).exp() / (s * (2.0 * std::f64::consts::PI).sqrt())
648}
649
650/// Assign functional group label to an IR frequency based on standard ranges.
651fn assign_ir_peak(frequency_cm1: f64) -> String {
652    if frequency_cm1 > 3500.0 {
653        "O-H stretch (broad)".to_string()
654    } else if frequency_cm1 > 3300.0 {
655        "N-H stretch".to_string()
656    } else if frequency_cm1 > 3000.0 {
657        "sp2 C-H stretch".to_string()
658    } else if frequency_cm1 > 2800.0 {
659        "sp3 C-H stretch".to_string()
660    } else if frequency_cm1 > 2100.0 && frequency_cm1 < 2300.0 {
661        "C≡N or C≡C stretch".to_string()
662    } else if frequency_cm1 > 1680.0 && frequency_cm1 < 1800.0 {
663        "C=O stretch".to_string()
664    } else if frequency_cm1 > 1550.0 && frequency_cm1 < 1680.0 {
665        "H-O-H bend / N-H bend".to_string()
666    } else if frequency_cm1 > 1600.0 && frequency_cm1 < 1680.0 {
667        "C=C stretch".to_string()
668    } else if frequency_cm1 > 1400.0 && frequency_cm1 < 1600.0 {
669        "aromatic C=C stretch".to_string()
670    } else if frequency_cm1 > 1300.0 && frequency_cm1 < 1400.0 {
671        "C-H bend".to_string()
672    } else if frequency_cm1 > 1000.0 && frequency_cm1 < 1300.0 {
673        "C-O stretch".to_string()
674    } else if frequency_cm1 > 600.0 && frequency_cm1 < 900.0 {
675        "C-H out-of-plane bend".to_string()
676    } else {
677        "skeletal mode".to_string()
678    }
679}
680
681/// Compute RRHO thermochemistry from vibrational frequencies.
682///
683/// Uses the rigid-rotor harmonic-oscillator approximation at the given temperature.
684fn compute_thermochemistry(frequencies_cm1: &[f64], temperature_k: f64) -> Thermochemistry {
685    // Constants
686    const CM1_TO_EV: f64 = 1.0 / 8065.544;
687    const EV_TO_KCAL: f64 = 23.0605;
688    const R_KCAL: f64 = 0.001987204; // kcal/(mol·K)
689    const R_CAL: f64 = 1.987204; // cal/(mol·K)
690    const KB_EV: f64 = 8.617333e-5; // eV/K
691
692    let kbt_ev = KB_EV * temperature_k;
693
694    let real_freqs: Vec<f64> = frequencies_cm1
695        .iter()
696        .filter(|&&f| f > 50.0)
697        .copied()
698        .collect();
699
700    // ZPVE = Σ (1/2)hν with 0.9 anharmonic scaling
701    let zpve_ev: f64 = real_freqs.iter().map(|&f| 0.5 * f * CM1_TO_EV * 0.9).sum();
702    let zpve_kcal = zpve_ev * EV_TO_KCAL;
703
704    // Thermal energy (vibrational contribution): Σ hν / (exp(hν/kT) - 1)
705    let thermal_e_ev: f64 = real_freqs
706        .iter()
707        .map(|&f| {
708            let hnu = f * CM1_TO_EV;
709            let x = hnu / kbt_ev;
710            if x > 100.0 {
711                0.0
712            } else {
713                hnu / (x.exp() - 1.0)
714            }
715        })
716        .sum();
717    let thermal_energy_kcal = thermal_e_ev * EV_TO_KCAL;
718
719    // Enthalpy correction: E_thermal + RT
720    let enthalpy_correction_kcal = thermal_energy_kcal + R_KCAL * temperature_k;
721
722    // Vibrational entropy: Σ [x/(exp(x)-1) - ln(1-exp(-x))]·R
723    let entropy_vib_cal: f64 = real_freqs
724        .iter()
725        .map(|&f| {
726            let hnu = f * CM1_TO_EV;
727            let x = hnu / kbt_ev;
728            if x > 100.0 {
729                0.0
730            } else {
731                let ex = x.exp();
732                R_CAL * (x / (ex - 1.0) - (1.0 - (-x).exp()).ln())
733            }
734        })
735        .sum();
736
737    // Gibbs correction: H - TS
738    let gibbs_correction_kcal = enthalpy_correction_kcal - temperature_k * entropy_vib_cal / 1000.0;
739
740    Thermochemistry {
741        temperature_k,
742        zpve_kcal,
743        thermal_energy_kcal,
744        enthalpy_correction_kcal,
745        entropy_vib_cal,
746        gibbs_correction_kcal,
747    }
748}
749
750/// Perform a complete vibrational analysis from elements and positions.
751///
752/// Computes the numerical Hessian, mass-weighted eigenvalue problem,
753/// vibrational frequencies, and IR intensities.
754pub fn compute_vibrational_analysis(
755    elements: &[u8],
756    positions: &[[f64; 3]],
757    method: HessianMethod,
758    step_size: Option<f64>,
759) -> Result<VibrationalAnalysis, String> {
760    if method == HessianMethod::Uff {
761        return Err(
762            "UFF requires SMILES; use compute_vibrational_analysis_uff instead".to_string(),
763        );
764    }
765
766    let n_atoms = elements.len();
767    let delta = step_size.unwrap_or(0.005);
768
769    if n_atoms < 2 {
770        return Err("Need at least 2 atoms for vibrational analysis".to_string());
771    }
772
773    let hessian = compute_numerical_hessian(elements, positions, method, Some(delta))?;
774    build_vibrational_analysis_from_hessian(elements, positions, &hessian, method, delta, None)
775}
776
777/// Perform vibrational analysis using UFF analytical Hessian (gradient-difference method).
778///
779/// This avoids the expensive O(9N²) energy evaluations of the standard numerical Hessian
780/// by using O(6N) gradient evaluations with UFF's analytical gradients.
781pub fn compute_vibrational_analysis_uff(
782    smiles: &str,
783    elements: &[u8],
784    positions: &[[f64; 3]],
785    step_size: Option<f64>,
786) -> Result<VibrationalAnalysis, String> {
787    let n_atoms = elements.len();
788    let delta = step_size.unwrap_or(0.005);
789
790    if n_atoms < 2 {
791        return Err("Need at least 2 atoms for vibrational analysis".to_string());
792    }
793
794    let relaxed_positions = relax_uff_geometry(smiles, positions)?;
795    let coords_flat: Vec<f64> = relaxed_positions
796        .iter()
797        .flat_map(|position| position.iter().copied())
798        .collect();
799    let hessian =
800        super::hessian::compute_uff_analytical_hessian(smiles, &coords_flat, Some(delta))?;
801    let uff_charges = crate::compute_charges(smiles)
802        .ok()
803        .map(|result| result.charges);
804    let mut analysis = build_vibrational_analysis_from_hessian(
805        elements,
806        &relaxed_positions,
807        &hessian,
808        HessianMethod::Uff,
809        delta,
810        uff_charges.as_deref(),
811    )?;
812
813    if has_significant_imaginary_vibrational_mode(&analysis) {
814        let numerical_hessian = compute_uff_numerical_hessian(smiles, &coords_flat, delta)?;
815        analysis = build_vibrational_analysis_from_hessian(
816            elements,
817            &relaxed_positions,
818            &numerical_hessian,
819            HessianMethod::Uff,
820            delta,
821            uff_charges.as_deref(),
822        )?;
823        analysis.notes.push(
824            "UFF analytical Hessian retained significant imaginary vibrational modes, so the analysis fell back to a numerical energy Hessian for stability.".to_string(),
825        );
826    }
827
828    analysis
829        .notes
830        .push("UFF coordinates were pre-relaxed with an internal gradient descent before the Hessian was built.".to_string());
831    Ok(analysis)
832}
833
834/// Build vibrational analysis from a pre-computed Hessian matrix.
835fn build_vibrational_analysis_from_hessian(
836    elements: &[u8],
837    positions: &[[f64; 3]],
838    hessian: &DMatrix<f64>,
839    method: HessianMethod,
840    delta: f64,
841    uff_charges: Option<&[f64]>,
842) -> Result<VibrationalAnalysis, String> {
843    let n_atoms = elements.len();
844    let n3 = 3 * n_atoms;
845
846    // 2. Build mass-weighted Hessian: H'_{ij} = H_{ij} / sqrt(m_i * m_j)
847    let masses: Vec<f64> = elements.iter().map(|&z| atomic_mass(z)).collect();
848    let mut mw_hessian = DMatrix::zeros(n3, n3);
849    for i in 0..n3 {
850        let mi = masses[i / 3];
851        for j in 0..n3 {
852            let mj = masses[j / 3];
853            mw_hessian[(i, j)] = hessian[(i, j)] / (mi * mj).sqrt();
854        }
855    }
856
857    let rigid_basis = build_rigid_body_basis(&masses, positions);
858    let projected_hessian = project_rigid_body_modes(&mw_hessian, &rigid_basis);
859
860    // 3. Diagonalize mass-weighted Hessian
861    let eigen = projected_hessian.symmetric_eigen();
862
863    // 4. Sort eigenvalues and eigenvectors by eigenvalue
864    let mut indices: Vec<usize> = (0..n3).collect();
865    indices.sort_by(|&a, &b| {
866        eigen.eigenvalues[a]
867            .partial_cmp(&eigen.eigenvalues[b])
868            .unwrap_or(std::cmp::Ordering::Equal)
869    });
870
871    // 5. Convert eigenvalues to frequencies (cm⁻¹)
872    let frequency_factor = hessian_frequency_factor(method);
873
874    let mut sorted_eigenvalues = Vec::with_capacity(n3);
875    let mut sorted_modes = DMatrix::zeros(n3, n3);
876    for (new_idx, &old_idx) in indices.iter().enumerate() {
877        sorted_eigenvalues.push(eigen.eigenvalues[old_idx]);
878        for i in 0..n3 {
879            sorted_modes[(i, new_idx)] = eigen.eigenvectors[(i, old_idx)];
880        }
881    }
882
883    // Convert eigenvalues to frequencies
884    let frequencies: Vec<f64> = sorted_eigenvalues
885        .iter()
886        .map(|&ev| {
887            if ev >= 0.0 {
888                // ν = (1/2πc) * sqrt(eigenvalue * factor)
889                (ev * frequency_factor).sqrt() * INV_2PI_C
890            } else {
891                // Imaginary frequency (negative eigenvalue)
892                -((-ev) * frequency_factor).sqrt() * INV_2PI_C
893            }
894        })
895        .collect();
896
897    // 6. Compute IR intensities from dipole derivatives
898    let ir_intensities = compute_ir_intensities(
899        elements,
900        positions,
901        &sorted_modes,
902        method,
903        delta * 2.0, // slightly larger step for dipole derivatives
904        uff_charges,
905    )?;
906
907    // 7. Build mode list
908    let mut modes = Vec::with_capacity(n3);
909    let mut n_real = 0;
910    let mut zpve = 0.0;
911
912    for k in 0..n3 {
913        let freq = frequencies[k];
914        let displacement: Vec<f64> = (0..n3).map(|i| sorted_modes[(i, k)]).collect();
915        let is_real = rigid_body_overlap(&displacement, &rigid_basis) < 0.5;
916
917        if is_real {
918            n_real += 1;
919            if freq > 0.0 {
920                // ZPVE contribution: (1/2)hν in eV
921                zpve += 0.5 * freq / EV_TO_CM1;
922            }
923        }
924
925        modes.push(VibrationalMode {
926            frequency_cm1: freq,
927            ir_intensity: ir_intensities.get(k).copied().unwrap_or(0.0),
928            displacement,
929            is_real,
930        });
931    }
932
933    let method_name = match method {
934        HessianMethod::Eht => "EHT",
935        HessianMethod::Pm3 => "PM3",
936        HessianMethod::Xtb => "xTB",
937        HessianMethod::Uff => "UFF",
938    };
939
940    // Compute thermochemistry at 298.15 K
941    let thermochemistry = Some(compute_thermochemistry(&frequencies, 298.15));
942
943    Ok(VibrationalAnalysis {
944        n_atoms,
945        elements: elements.to_vec(),
946        modes,
947        n_real_modes: n_real,
948        zpve_ev: zpve,
949        thermochemistry,
950        method: method_name.to_string(),
951        notes: vec![
952            format!(
953                "Numerical Hessian computed with {} using central finite differences (δ = {} Å).",
954                method_name, delta
955            ),
956            match method {
957                HessianMethod::Uff if uff_charges.is_some() => {
958                    "IR intensities derived from numerical dipole derivatives using a Gasteiger-charge dipole approximation.".to_string()
959                }
960                HessianMethod::Uff => {
961                    "IR intensities derived from a neutral-charge dipole approximation because Gasteiger charges were unavailable.".to_string()
962                }
963                _ => "IR intensities derived from numerical dipole derivatives along normal coordinates.".to_string(),
964            },
965            format!(
966                "Rigid-body translation/rotation modes were projected out before diagonalization ({} basis vectors removed).",
967                rigid_basis.ncols()
968            ),
969            match method {
970                HessianMethod::Uff => {
971                    "UFF Hessian values were converted from kcal/mol·Å⁻² to the shared frequency scale before computing cm⁻¹ bands.".to_string()
972                }
973                _ => "Semiempirical Hessian values were interpreted on the shared eV·Å⁻² frequency scale.".to_string(),
974            },
975        ],
976    })
977}
978
979/// Broadening function type for IR spectrum generation.
980#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
981pub enum BroadeningType {
982    /// Lorentzian line shape (default).
983    Lorentzian,
984    /// Gaussian line shape.
985    Gaussian,
986}
987
988/// Generate a broadened IR spectrum from vibrational analysis.
989///
990/// `gamma`: broadening width in cm⁻¹ (typically 10–30).
991/// `wn_min`, `wn_max`: wavenumber range in cm⁻¹.
992/// `n_points`: number of grid points.
993pub fn compute_ir_spectrum(
994    analysis: &VibrationalAnalysis,
995    gamma: f64,
996    wn_min: f64,
997    wn_max: f64,
998    n_points: usize,
999) -> IrSpectrum {
1000    compute_ir_spectrum_with_broadening(
1001        analysis,
1002        gamma,
1003        wn_min,
1004        wn_max,
1005        n_points,
1006        BroadeningType::Lorentzian,
1007    )
1008}
1009
1010/// Generate a broadened IR spectrum with selectable broadening function.
1011pub fn compute_ir_spectrum_with_broadening(
1012    analysis: &VibrationalAnalysis,
1013    gamma: f64,
1014    wn_min: f64,
1015    wn_max: f64,
1016    n_points: usize,
1017    broadening: BroadeningType,
1018) -> IrSpectrum {
1019    let n_points = n_points.max(2);
1020    let step = (wn_max - wn_min) / (n_points as f64 - 1.0);
1021    let wavenumbers: Vec<f64> = (0..n_points).map(|i| wn_min + step * i as f64).collect();
1022    let mut intensities = vec![0.0; n_points];
1023
1024    let active_modes: Vec<(usize, &VibrationalMode)> = analysis
1025        .modes
1026        .iter()
1027        .enumerate()
1028        .filter(|(_, mode)| mode.is_real && mode.frequency_cm1 > 0.0)
1029        .collect();
1030    let active_frequencies: Vec<f64> = active_modes
1031        .iter()
1032        .map(|(_, mode)| mode.frequency_cm1)
1033        .collect();
1034    let active_mode_intensities: Vec<f64> = active_modes
1035        .iter()
1036        .map(|(_, mode)| mode.ir_intensity)
1037        .collect();
1038    let active_mode_displacements: Vec<Vec<[f64; 3]>> = active_modes
1039        .iter()
1040        .map(|(_, mode)| {
1041            mode.displacement
1042                .chunks(3)
1043                .map(|chunk| [chunk[0], chunk[1], chunk[2]])
1044                .collect()
1045        })
1046        .collect();
1047    let assigned_labels = if analysis.elements.is_empty() {
1048        vec![None; active_modes.len()]
1049    } else {
1050        let assignment_result = super::peak_assignment::assign_peaks(
1051            &active_frequencies,
1052            &active_mode_intensities,
1053            &analysis.elements,
1054            Some(&active_mode_displacements),
1055        );
1056        let mut labels = vec![None; active_modes.len()];
1057        let mut assigned_idx = 0usize;
1058
1059        for (mode_idx, frequency) in active_frequencies.iter().enumerate() {
1060            if *frequency < 400.0 {
1061                continue;
1062            }
1063            if let Some(assignment) = assignment_result.assignments.get(assigned_idx) {
1064                if let Some(best) = assignment.assignments.iter().max_by(|left, right| {
1065                    left.match_quality
1066                        .partial_cmp(&right.match_quality)
1067                        .unwrap_or(std::cmp::Ordering::Equal)
1068                }) {
1069                    labels[mode_idx] = Some(if best.group.contains(&best.vibration_type) {
1070                        best.group.clone()
1071                    } else {
1072                        format!("{} {}", best.group, best.vibration_type)
1073                    });
1074                }
1075            }
1076            assigned_idx += 1;
1077        }
1078
1079        labels
1080    };
1081
1082    let mut peaks = Vec::new();
1083
1084    for (active_idx, (mode_idx, mode)) in active_modes.iter().enumerate() {
1085        peaks.push(IrPeak {
1086            frequency_cm1: mode.frequency_cm1,
1087            ir_intensity: mode.ir_intensity,
1088            mode_index: *mode_idx,
1089            assignment: assigned_labels
1090                .get(active_idx)
1091                .and_then(|label| label.clone())
1092                .unwrap_or_else(|| assign_ir_peak(mode.frequency_cm1)),
1093        });
1094
1095        for (idx, &wn) in wavenumbers.iter().enumerate() {
1096            intensities[idx] += mode.ir_intensity
1097                * match broadening {
1098                    BroadeningType::Lorentzian => lorentzian(wn, mode.frequency_cm1, gamma),
1099                    BroadeningType::Gaussian => gaussian(wn, mode.frequency_cm1, gamma),
1100                };
1101        }
1102    }
1103
1104    peaks.sort_by(|a, b| {
1105        b.ir_intensity
1106            .partial_cmp(&a.ir_intensity)
1107            .unwrap_or(std::cmp::Ordering::Equal)
1108    });
1109
1110    IrSpectrum {
1111        wavenumbers,
1112        intensities,
1113        peaks,
1114        gamma,
1115        notes: vec![
1116            format!(
1117                "IR spectrum generated with Lorentzian broadening (γ = {} cm⁻¹).",
1118                gamma
1119            ),
1120            format!("Vibrational analysis method: {}.", analysis.method),
1121        ],
1122    }
1123}
1124
1125#[cfg(test)]
1126mod tests {
1127    use super::*;
1128
1129    #[test]
1130    fn test_lorentzian_peak_at_center() {
1131        let val = lorentzian(100.0, 100.0, 10.0);
1132        // At center: L = 1/(π·γ)
1133        let expected = 1.0 / (std::f64::consts::PI * 10.0);
1134        assert!((val - expected).abs() < 1e-10);
1135    }
1136
1137    #[test]
1138    fn test_lorentzian_symmetry() {
1139        let left = lorentzian(90.0, 100.0, 10.0);
1140        let right = lorentzian(110.0, 100.0, 10.0);
1141        assert!((left - right).abs() < 1e-10);
1142    }
1143
1144    #[test]
1145    fn test_atomic_masses_known_elements() {
1146        assert!((atomic_mass(1) - 1.00794).abs() < 0.001);
1147        assert!((atomic_mass(6) - 12.011).abs() < 0.01);
1148        assert!((atomic_mass(8) - 15.999).abs() < 0.01);
1149    }
1150
1151    #[test]
1152    fn test_h2_vibrational_analysis() {
1153        let elements = [1u8, 1];
1154        let positions = [[0.0, 0.0, 0.0], [0.74, 0.0, 0.0]];
1155        let analysis =
1156            compute_vibrational_analysis(&elements, &positions, HessianMethod::Xtb, Some(0.005))
1157                .unwrap();
1158
1159        assert_eq!(analysis.n_atoms, 2);
1160        assert_eq!(analysis.modes.len(), 6); // 3N = 6
1161
1162        // H₂ is linear: 5 zero + 1 stretch
1163        // Should have at least 1 real mode (H-H stretch)
1164        assert!(
1165            analysis.n_real_modes >= 1,
1166            "H₂ should have at least 1 real vibrational mode, got {}",
1167            analysis.n_real_modes
1168        );
1169
1170        // The stretch frequency should be positive
1171        let real_modes: Vec<&VibrationalMode> =
1172            analysis.modes.iter().filter(|m| m.is_real).collect();
1173        assert!(!real_modes.is_empty(), "Should have at least one real mode");
1174
1175        // Check that ZPVE is positive for real modes
1176        if analysis.n_real_modes > 0 {
1177            assert!(analysis.zpve_ev > 0.0, "ZPVE should be positive");
1178        }
1179    }
1180
1181    #[test]
1182    fn test_water_vibrational_analysis() {
1183        let elements = [8u8, 1, 1];
1184        let positions = [[0.0, 0.0, 0.0], [0.757, 0.586, 0.0], [-0.757, 0.586, 0.0]];
1185        let analysis =
1186            compute_vibrational_analysis(&elements, &positions, HessianMethod::Xtb, Some(0.005))
1187                .unwrap();
1188
1189        assert_eq!(analysis.n_atoms, 3);
1190        assert_eq!(analysis.modes.len(), 9); // 3N = 9
1191
1192        // Water is nonlinear: 6 zero + 3 real modes (sym stretch, asym stretch, bend)
1193        // Due to numerical noise some translational modes may appear as real
1194        // so we just check that we have at least 3 real modes
1195        assert!(
1196            analysis.n_real_modes >= 2,
1197            "Water should have at least 2-3 real modes, got {}",
1198            analysis.n_real_modes
1199        );
1200    }
1201
1202    #[test]
1203    fn test_ir_spectrum_generation() {
1204        let elements = [8u8, 1, 1];
1205        let positions = [[0.0, 0.0, 0.0], [0.757, 0.586, 0.0], [-0.757, 0.586, 0.0]];
1206        let analysis =
1207            compute_vibrational_analysis(&elements, &positions, HessianMethod::Xtb, Some(0.005))
1208                .unwrap();
1209
1210        let spectrum = compute_ir_spectrum(&analysis, 20.0, 400.0, 4000.0, 500);
1211
1212        assert_eq!(spectrum.wavenumbers.len(), 500);
1213        assert_eq!(spectrum.intensities.len(), 500);
1214        assert!(!spectrum.peaks.is_empty(), "Should have IR peaks");
1215        assert!(
1216            spectrum.intensities.iter().any(|&i| i > 0.0),
1217            "Spectrum should have non-zero intensity"
1218        );
1219    }
1220}