Skip to main content

sci_form/
dynamics.rs

1use rand::rngs::StdRng;
2use rand::{Rng, SeedableRng};
3use serde::{Deserialize, Serialize};
4
5const AMU_ANGFS2_TO_KCAL_MOL: f64 = 2_390.057_361_533_49;
6const R_GAS_KCAL_MOLK: f64 = 0.001_987_204_258_640_83;
7const EV_TO_KCAL_MOL: f64 = 23.060_5;
8const HARTREE_TO_KCAL_MOL: f64 = 627.509_474_063_1;
9
10/// Force backend for molecular dynamics.
11#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
12pub enum MdBackend {
13    /// UFF force field (fast, approximate).
14    Uff,
15    /// PM3 semi-empirical (slower, more accurate).
16    Pm3,
17    /// GFN-xTB tight-binding (moderate speed, good for metals).
18    Xtb,
19}
20
21/// Energy backend for NEB path calculations.
22///
23/// Supports all methods that can provide energy and gradients (analytical
24/// or numerical) for NEB image relaxation.
25#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
26#[serde(rename_all = "snake_case")]
27pub enum NebBackend {
28    /// UFF force field — fast exploratory paths.
29    Uff,
30    /// MMFF94 force field — better organic coverage than UFF.
31    Mmff94,
32    /// PM3 semi-empirical — analytical gradients.
33    Pm3,
34    /// GFN0-xTB tight-binding — analytical gradients.
35    Xtb,
36    /// GFN1-xTB — numerical gradients (energy-only + finite differences).
37    Gfn1,
38    /// GFN2-xTB — numerical gradients (energy-only + finite differences).
39    Gfn2,
40    /// HF-3c minimal-basis Hartree-Fock — numerical gradients (expensive).
41    Hf3c,
42}
43
44impl NebBackend {
45    /// Parse a method string into a `NebBackend`.
46    pub fn from_method(s: &str) -> Result<Self, String> {
47        match s.trim().to_ascii_lowercase().as_str() {
48            "uff" => Ok(Self::Uff),
49            "mmff94" | "mmff" => Ok(Self::Mmff94),
50            "pm3" => Ok(Self::Pm3),
51            "xtb" | "gfn0" | "gfn0-xtb" | "gfn0_xtb" => Ok(Self::Xtb),
52            "gfn1" | "gfn1-xtb" | "gfn1_xtb" => Ok(Self::Gfn1),
53            "gfn2" | "gfn2-xtb" | "gfn2_xtb" => Ok(Self::Gfn2),
54            "hf3c" | "hf-3c" => Ok(Self::Hf3c),
55            other => Err(format!(
56                "Unknown NEB backend '{}'. Expected: uff, mmff94, pm3, xtb, gfn1, gfn2, hf3c",
57                other
58            )),
59        }
60    }
61
62    /// method label for human display.
63    pub fn as_str(self) -> &'static str {
64        match self {
65            Self::Uff => "uff",
66            Self::Mmff94 => "mmff94",
67            Self::Pm3 => "pm3",
68            Self::Xtb => "xtb",
69            Self::Gfn1 => "gfn1",
70            Self::Gfn2 => "gfn2",
71            Self::Hf3c => "hf3c",
72        }
73    }
74}
75
76/// One trajectory frame for molecular-dynamics sampling.
77#[derive(Debug, Clone, Serialize, Deserialize)]
78pub struct MdFrame {
79    /// Integration step index.
80    pub step: usize,
81    /// Elapsed simulation time in femtoseconds.
82    pub time_fs: f64,
83    /// Flat xyz coordinates in angstroms.
84    pub coords: Vec<f64>,
85    /// Potential energy from UFF in kcal/mol.
86    pub potential_energy_kcal_mol: f64,
87    /// Kinetic energy in kcal/mol.
88    pub kinetic_energy_kcal_mol: f64,
89    /// Instantaneous temperature estimate in K.
90    pub temperature_k: f64,
91}
92
93/// Full trajectory output for exploratory molecular dynamics.
94#[derive(Debug, Clone, Serialize, Deserialize)]
95pub struct MdTrajectory {
96    /// Stored MD frames.
97    pub frames: Vec<MdFrame>,
98    /// Timestep in femtoseconds.
99    pub dt_fs: f64,
100    /// Notes and caveats for interpretation.
101    pub notes: Vec<String>,
102    /// Energy conservation drift: (E_total_final - E_total_initial) / |E_total_initial| * 100%.
103    /// Only meaningful for NVE (no thermostat) simulations.
104    pub energy_drift_percent: Option<f64>,
105}
106
107/// One image (node) on a simplified NEB path.
108#[derive(Debug, Clone, Serialize, Deserialize)]
109pub struct NebImage {
110    /// Image index from reactant (0) to product (n-1).
111    pub index: usize,
112    /// Flat xyz coordinates in angstroms.
113    pub coords: Vec<f64>,
114    /// UFF potential energy in kcal/mol.
115    pub potential_energy_kcal_mol: f64,
116}
117
118/// Simplified NEB output for low-cost pathway exploration.
119#[derive(Debug, Clone, Serialize, Deserialize)]
120pub struct NebPathResult {
121    /// Ordered path images from reactant to product.
122    pub images: Vec<NebImage>,
123    /// Notes and caveats.
124    pub notes: Vec<String>,
125}
126
127pub fn atomic_mass_amu(z: u8) -> f64 {
128    match z {
129        1 => 1.008,
130        5 => 10.81,
131        6 => 12.011,
132        7 => 14.007,
133        8 => 15.999,
134        9 => 18.998,
135        14 => 28.085,
136        15 => 30.974,
137        16 => 32.06,
138        17 => 35.45,
139        35 => 79.904,
140        53 => 126.904,
141        26 => 55.845,
142        46 => 106.42,
143        78 => 195.084,
144        _ => 12.0,
145    }
146}
147
148fn sample_standard_normal(rng: &mut StdRng) -> f64 {
149    let u1 = (1.0 - rng.gen::<f64>()).max(1e-12);
150    let u2 = rng.gen::<f64>();
151    (-2.0 * u1.ln()).sqrt() * (2.0 * std::f64::consts::PI * u2).cos()
152}
153
154fn kinetic_energy_and_temperature(
155    velocities: &[f64],
156    masses_amu: &[f64],
157    n_atoms: usize,
158) -> (f64, f64) {
159    let mut ke = 0.0;
160    for i in 0..n_atoms {
161        let vx = velocities[3 * i];
162        let vy = velocities[3 * i + 1];
163        let vz = velocities[3 * i + 2];
164        let v2 = vx * vx + vy * vy + vz * vz;
165        ke += 0.5 * masses_amu[i] * v2 * AMU_ANGFS2_TO_KCAL_MOL;
166    }
167    let dof = (3 * n_atoms).saturating_sub(6).max(1) as f64;
168    let t = 2.0 * ke / (dof * R_GAS_KCAL_MOLK);
169    (ke, t)
170}
171
172/// Run short exploratory molecular dynamics using Velocity Verlet and optional Berendsen NVT.
173pub fn simulate_velocity_verlet_uff(
174    smiles: &str,
175    coords: &[f64],
176    n_steps: usize,
177    dt_fs: f64,
178    seed: u64,
179    target_temp_and_tau: Option<(f64, f64)>,
180) -> Result<MdTrajectory, String> {
181    if n_steps == 0 {
182        return Err("n_steps must be > 0".to_string());
183    }
184    if dt_fs <= 0.0 {
185        return Err("dt_fs must be > 0".to_string());
186    }
187
188    let mol = crate::graph::Molecule::from_smiles(smiles)?;
189    let n_atoms = mol.graph.node_count();
190    if coords.len() != n_atoms * 3 {
191        return Err(format!(
192            "coords length {} != 3 * atoms {}",
193            coords.len(),
194            n_atoms
195        ));
196    }
197
198    let masses_amu: Vec<f64> = (0..n_atoms)
199        .map(petgraph::graph::NodeIndex::new)
200        .map(|idx| atomic_mass_amu(mol.graph[idx].element))
201        .collect();
202
203    let ff = crate::forcefield::builder::build_uff_force_field(&mol);
204    let mut x = coords.to_vec();
205    let mut grad = vec![0.0f64; n_atoms * 3];
206    let mut potential = ff.compute_system_energy_and_gradients(&x, &mut grad);
207
208    let mut rng = StdRng::seed_from_u64(seed);
209    let mut v = vec![0.0f64; n_atoms * 3];
210
211    if let Some((target_temp_k, _tau_fs)) = target_temp_and_tau {
212        for i in 0..n_atoms {
213            let sigma = ((R_GAS_KCAL_MOLK * target_temp_k)
214                / (masses_amu[i] * AMU_ANGFS2_TO_KCAL_MOL))
215                .sqrt();
216            v[3 * i] = sigma * sample_standard_normal(&mut rng);
217            v[3 * i + 1] = sigma * sample_standard_normal(&mut rng);
218            v[3 * i + 2] = sigma * sample_standard_normal(&mut rng);
219        }
220    }
221
222    let (ke0, t0) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
223    let mut frames = vec![MdFrame {
224        step: 0,
225        time_fs: 0.0,
226        coords: x.clone(),
227        potential_energy_kcal_mol: potential,
228        kinetic_energy_kcal_mol: ke0,
229        temperature_k: t0,
230    }];
231
232    for step in 1..=n_steps {
233        let mut a = vec![0.0f64; n_atoms * 3];
234        for i in 0..n_atoms {
235            let m = masses_amu[i];
236            a[3 * i] = -grad[3 * i] / (m * AMU_ANGFS2_TO_KCAL_MOL);
237            a[3 * i + 1] = -grad[3 * i + 1] / (m * AMU_ANGFS2_TO_KCAL_MOL);
238            a[3 * i + 2] = -grad[3 * i + 2] / (m * AMU_ANGFS2_TO_KCAL_MOL);
239        }
240
241        for i in 0..(n_atoms * 3) {
242            x[i] += v[i] * dt_fs + 0.5 * a[i] * dt_fs * dt_fs;
243        }
244
245        potential = ff.compute_system_energy_and_gradients(&x, &mut grad);
246
247        let mut a_new = vec![0.0f64; n_atoms * 3];
248        for i in 0..n_atoms {
249            let m = masses_amu[i];
250            a_new[3 * i] = -grad[3 * i] / (m * AMU_ANGFS2_TO_KCAL_MOL);
251            a_new[3 * i + 1] = -grad[3 * i + 1] / (m * AMU_ANGFS2_TO_KCAL_MOL);
252            a_new[3 * i + 2] = -grad[3 * i + 2] / (m * AMU_ANGFS2_TO_KCAL_MOL);
253        }
254
255        for i in 0..(n_atoms * 3) {
256            v[i] += 0.5 * (a[i] + a_new[i]) * dt_fs;
257        }
258
259        let (mut ke, mut temp_k) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
260        if let Some((target_temp_k, tau_fs)) = target_temp_and_tau {
261            let tau = tau_fs.max(1e-6);
262            let lambda = (1.0 + (dt_fs / tau) * (target_temp_k / temp_k.max(1e-6) - 1.0)).sqrt();
263            let lambda = lambda.clamp(0.5, 2.0);
264            for vi in &mut v {
265                *vi *= lambda;
266            }
267            let kt = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
268            ke = kt.0;
269            temp_k = kt.1;
270        }
271
272        if !x.iter().all(|v| v.is_finite()) || !potential.is_finite() || !ke.is_finite() {
273            return Err(format!(
274                "MD diverged at step {} (non-finite coordinates/energy)",
275                step
276            ));
277        }
278
279        frames.push(MdFrame {
280            step,
281            time_fs: step as f64 * dt_fs,
282            coords: x.clone(),
283            potential_energy_kcal_mol: potential,
284            kinetic_energy_kcal_mol: ke,
285            temperature_k: temp_k,
286        });
287    }
288
289    let mut notes = vec![
290        "Velocity Verlet integration over UFF force-field gradients for short exploratory trajectories."
291            .to_string(),
292    ];
293    if target_temp_and_tau.is_some() {
294        notes.push(
295            "Berendsen thermostat rescaling applied for approximate constant-temperature sampling."
296                .to_string(),
297        );
298    } else {
299        notes.push(
300            "No thermostat applied (NVE-like propagation with current numerical approximations)."
301                .to_string(),
302        );
303    }
304
305    Ok(MdTrajectory {
306        energy_drift_percent: if target_temp_and_tau.is_none() {
307            let e0 = frames[0].potential_energy_kcal_mol + frames[0].kinetic_energy_kcal_mol;
308            let ef = frames
309                .last()
310                .map(|f| f.potential_energy_kcal_mol + f.kinetic_energy_kcal_mol)
311                .unwrap_or(e0);
312            if e0.abs() > 1e-10 {
313                Some((ef - e0) / e0.abs() * 100.0)
314            } else {
315                None
316            }
317        } else {
318            None
319        },
320        frames,
321        dt_fs,
322        notes,
323    })
324}
325
326/// Build a simplified NEB-like path using linear interpolation and spring-coupled relaxation.
327pub fn compute_simplified_neb_path(
328    smiles: &str,
329    start_coords: &[f64],
330    end_coords: &[f64],
331    n_images: usize,
332    n_iter: usize,
333    spring_k: f64,
334    step_size: f64,
335) -> Result<NebPathResult, String> {
336    if n_images < 2 {
337        return Err("n_images must be >= 2".to_string());
338    }
339    if n_iter == 0 {
340        return Err("n_iter must be > 0".to_string());
341    }
342    if step_size <= 0.0 {
343        return Err("step_size must be > 0".to_string());
344    }
345
346    let mol = crate::graph::Molecule::from_smiles(smiles)?;
347    let n_atoms = mol.graph.node_count();
348    let n_xyz = n_atoms * 3;
349    if start_coords.len() != n_xyz || end_coords.len() != n_xyz {
350        return Err(format!(
351            "start/end coords must each have length {} (3 * n_atoms)",
352            n_xyz
353        ));
354    }
355
356    let ff = crate::forcefield::builder::build_uff_force_field(&mol);
357
358    let mut images = vec![vec![0.0f64; n_xyz]; n_images];
359    for (img_idx, img) in images.iter_mut().enumerate() {
360        let t = img_idx as f64 / (n_images - 1) as f64;
361        for k in 0..n_xyz {
362            img[k] = (1.0 - t) * start_coords[k] + t * end_coords[k];
363        }
364    }
365
366    for _ in 0..n_iter {
367        let prev = images.clone();
368
369        #[cfg(feature = "parallel")]
370        {
371            use rayon::prelude::*;
372            let updated: Vec<(usize, Vec<f64>)> = (1..(n_images - 1))
373                .into_par_iter()
374                .map(|i| {
375                    let ff_local = crate::forcefield::builder::build_uff_force_field(&mol);
376                    let mut grad = vec![0.0f64; n_xyz];
377                    let _ = ff_local.compute_system_energy_and_gradients(&prev[i], &mut grad);
378                    let mut new_img = prev[i].clone();
379                    for k in 0..n_xyz {
380                        let spring_force =
381                            spring_k * (prev[i + 1][k] - 2.0 * prev[i][k] + prev[i - 1][k]);
382                        let total_force = -grad[k] + spring_force;
383                        new_img[k] = prev[i][k] + step_size * total_force;
384                    }
385                    (i, new_img)
386                })
387                .collect();
388            for (i, img) in updated {
389                images[i] = img;
390            }
391        }
392
393        #[cfg(not(feature = "parallel"))]
394        {
395            for i in 1..(n_images - 1) {
396                let mut grad = vec![0.0f64; n_xyz];
397                let _ = ff.compute_system_energy_and_gradients(&prev[i], &mut grad);
398                for k in 0..n_xyz {
399                    let spring_force =
400                        spring_k * (prev[i + 1][k] - 2.0 * prev[i][k] + prev[i - 1][k]);
401                    let total_force = -grad[k] + spring_force;
402                    images[i][k] = prev[i][k] + step_size * total_force;
403                }
404            }
405        }
406    }
407
408    let mut out_images = Vec::with_capacity(n_images);
409    for (i, coords) in images.into_iter().enumerate() {
410        let mut grad = vec![0.0f64; n_xyz];
411        let e = ff.compute_system_energy_and_gradients(&coords, &mut grad);
412        out_images.push(NebImage {
413            index: i,
414            coords,
415            potential_energy_kcal_mol: e,
416        });
417    }
418
419    Ok(NebPathResult {
420        images: out_images,
421        notes: vec![
422            "Simplified NEB: linear interpolation + spring-coupled UFF gradient relaxation on internal images."
423                .to_string(),
424            "This is a low-cost exploratory path tool and not a full climbing-image / tangent-projected NEB implementation."
425                .to_string(),
426        ],
427    })
428}
429
430// ─── Configurable NEB backend dispatch ──────────────────────────────────────
431
432/// Compute energy (kcal/mol) for a single point using a backend that only
433/// exposes energy (no analytical gradients). Used by the numerical gradient
434/// fallback for GFN1, GFN2, and HF-3c.
435fn neb_point_energy_kcal(
436    backend: NebBackend,
437    elements: &[u8],
438    coords: &[f64],
439) -> Result<f64, String> {
440    let positions: Vec<[f64; 3]> = coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
441    match backend {
442        NebBackend::Gfn1 => {
443            let r = crate::xtb::gfn1::solve_gfn1(elements, &positions)?;
444            Ok(r.total_energy * EV_TO_KCAL_MOL)
445        }
446        NebBackend::Gfn2 => {
447            let r = crate::xtb::gfn2::solve_gfn2(elements, &positions)?;
448            Ok(r.total_energy * EV_TO_KCAL_MOL)
449        }
450        NebBackend::Hf3c => {
451            let r = crate::hf::solve_hf3c(elements, &positions, &crate::hf::HfConfig::default())?;
452            Ok(r.energy * HARTREE_TO_KCAL_MOL)
453        }
454        _ => unreachable!("neb_point_energy_kcal only for energy-only backends"),
455    }
456}
457
458/// Compute energy and numerical gradient (central finite differences) for
459/// backends without analytical gradients.
460///
461/// When the `parallel` feature is enabled, displacement evaluations are
462/// distributed across the rayon thread pool for ~Nx speedup.
463fn neb_numerical_gradient(
464    backend: NebBackend,
465    elements: &[u8],
466    coords: &[f64],
467    grad: &mut [f64],
468) -> Result<f64, String> {
469    let delta = 1e-5; // Å
470    let e0 = neb_point_energy_kcal(backend, elements, coords)?;
471
472    #[cfg(feature = "parallel")]
473    {
474        use rayon::prelude::*;
475        let n = coords.len();
476        let results: Vec<Result<f64, String>> = (0..n)
477            .into_par_iter()
478            .map(|i| {
479                let mut displaced = coords.to_vec();
480                displaced[i] = coords[i] + delta;
481                let e_plus = neb_point_energy_kcal(backend, elements, &displaced)?;
482                displaced[i] = coords[i] - delta;
483                let e_minus = neb_point_energy_kcal(backend, elements, &displaced)?;
484                Ok((e_plus - e_minus) / (2.0 * delta))
485            })
486            .collect();
487        for (i, r) in results.into_iter().enumerate() {
488            grad[i] = r?;
489        }
490    }
491
492    #[cfg(not(feature = "parallel"))]
493    {
494        let mut displaced = coords.to_vec();
495        for i in 0..coords.len() {
496            displaced[i] = coords[i] + delta;
497            let e_plus = neb_point_energy_kcal(backend, elements, &displaced)?;
498            displaced[i] = coords[i] - delta;
499            let e_minus = neb_point_energy_kcal(backend, elements, &displaced)?;
500            displaced[i] = coords[i]; // restore
501            grad[i] = (e_plus - e_minus) / (2.0 * delta);
502        }
503    }
504
505    Ok(e0)
506}
507
508/// Compute energy (kcal/mol) and gradients (kcal/mol/Å) for a NEB image
509/// using the specified backend.
510pub fn neb_energy_and_gradient(
511    backend: NebBackend,
512    _smiles: &str,
513    elements: &[u8],
514    mol: &crate::graph::Molecule,
515    coords: &[f64],
516    grad: &mut [f64],
517) -> Result<f64, String> {
518    match backend {
519        NebBackend::Uff => {
520            let ff = crate::forcefield::builder::build_uff_force_field(mol);
521            Ok(ff.compute_system_energy_and_gradients(coords, grad))
522        }
523        NebBackend::Mmff94 => {
524            let bonds: Vec<(usize, usize, u8)> = mol
525                .graph
526                .edge_indices()
527                .map(|e| {
528                    let (a, b) = mol.graph.edge_endpoints(e).unwrap();
529                    let order = match mol.graph[e].order {
530                        crate::graph::BondOrder::Single => 1u8,
531                        crate::graph::BondOrder::Double => 2,
532                        crate::graph::BondOrder::Triple => 3,
533                        crate::graph::BondOrder::Aromatic => 2,
534                        crate::graph::BondOrder::Unknown => 1,
535                    };
536                    (a.index(), b.index(), order)
537                })
538                .collect();
539            let terms = crate::forcefield::mmff94::Mmff94Builder::build(elements, &bonds);
540            let (energy, g) =
541                crate::forcefield::mmff94::Mmff94Builder::total_energy(&terms, coords);
542            grad[..g.len()].copy_from_slice(&g);
543            Ok(energy)
544        }
545        NebBackend::Pm3 => {
546            let positions: Vec<[f64; 3]> = coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
547            let r = crate::pm3::gradients::compute_pm3_gradient(elements, &positions)?;
548            let energy_kcal = r.energy * EV_TO_KCAL_MOL;
549            for (a, g) in r.gradients.iter().enumerate() {
550                for d in 0..3 {
551                    grad[a * 3 + d] = g[d] * EV_TO_KCAL_MOL;
552                }
553            }
554            Ok(energy_kcal)
555        }
556        NebBackend::Xtb => {
557            let positions: Vec<[f64; 3]> = coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
558            let r = crate::xtb::gradients::compute_xtb_gradient(elements, &positions)?;
559            let energy_kcal = r.energy * EV_TO_KCAL_MOL;
560            for (a, g) in r.gradients.iter().enumerate() {
561                for d in 0..3 {
562                    grad[a * 3 + d] = g[d] * EV_TO_KCAL_MOL;
563                }
564            }
565            Ok(energy_kcal)
566        }
567        NebBackend::Gfn1 | NebBackend::Gfn2 | NebBackend::Hf3c => {
568            neb_numerical_gradient(backend, elements, coords, grad)
569        }
570    }
571}
572
573/// Build a simplified NEB path with a configurable energy backend.
574///
575/// This is the multi-method version of [`compute_simplified_neb_path`].
576/// Supply `method` as one of: `"uff"`, `"mmff94"`, `"pm3"`, `"xtb"`,
577/// `"gfn1"`, `"gfn2"`, `"hf3c"`.
578pub fn compute_simplified_neb_path_configurable(
579    smiles: &str,
580    start_coords: &[f64],
581    end_coords: &[f64],
582    n_images: usize,
583    n_iter: usize,
584    spring_k: f64,
585    step_size: f64,
586    method: &str,
587) -> Result<NebPathResult, String> {
588    let backend = NebBackend::from_method(method)?;
589    if n_images < 2 {
590        return Err("n_images must be >= 2".to_string());
591    }
592    if n_iter == 0 {
593        return Err("n_iter must be > 0".to_string());
594    }
595    if step_size <= 0.0 {
596        return Err("step_size must be > 0".to_string());
597    }
598
599    let mol = crate::graph::Molecule::from_smiles(smiles)?;
600    let n_atoms = mol.graph.node_count();
601    let n_xyz = n_atoms * 3;
602    if start_coords.len() != n_xyz || end_coords.len() != n_xyz {
603        return Err(format!(
604            "start/end coords must each have length {} (3 * n_atoms)",
605            n_xyz
606        ));
607    }
608
609    let elements: Vec<u8> = (0..n_atoms)
610        .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
611        .collect();
612
613    // Linear interpolation
614    let mut images = vec![vec![0.0f64; n_xyz]; n_images];
615    for (img_idx, img) in images.iter_mut().enumerate() {
616        let t = img_idx as f64 / (n_images - 1) as f64;
617        for k in 0..n_xyz {
618            img[k] = (1.0 - t) * start_coords[k] + t * end_coords[k];
619        }
620    }
621
622    // Spring-coupled NEB relaxation
623    for _ in 0..n_iter {
624        let prev = images.clone();
625        for i in 1..(n_images - 1) {
626            let mut grad = vec![0.0f64; n_xyz];
627            let _ = neb_energy_and_gradient(backend, smiles, &elements, &mol, &prev[i], &mut grad)?;
628            for k in 0..n_xyz {
629                let spring_force = spring_k * (prev[i + 1][k] - 2.0 * prev[i][k] + prev[i - 1][k]);
630                let total_force = -grad[k] + spring_force;
631                images[i][k] = prev[i][k] + step_size * total_force;
632            }
633        }
634    }
635
636    // Evaluate final energies
637    let mut out_images = Vec::with_capacity(n_images);
638    for (i, coords) in images.into_iter().enumerate() {
639        let mut grad = vec![0.0f64; n_xyz];
640        let e = neb_energy_and_gradient(backend, smiles, &elements, &mol, &coords, &mut grad)?;
641        out_images.push(NebImage {
642            index: i,
643            coords,
644            potential_energy_kcal_mol: e,
645        });
646    }
647
648    Ok(NebPathResult {
649        images: out_images,
650        notes: vec![
651            format!(
652                "Simplified NEB ({}) with spring-coupled gradient relaxation on {} internal images.",
653                backend.as_str(),
654                n_images.saturating_sub(2),
655            ),
656            "Low-cost exploratory path; not a full climbing-image / tangent-projected NEB."
657                .to_string(),
658        ],
659    })
660}
661
662/// Compute energy-only for any NEB backend (used for single-point comparisons).
663///
664/// Returns energy in kcal/mol.
665pub fn neb_backend_energy_kcal(method: &str, smiles: &str, coords: &[f64]) -> Result<f64, String> {
666    let backend = NebBackend::from_method(method)?;
667    let mol = crate::graph::Molecule::from_smiles(smiles)?;
668    let n_atoms = mol.graph.node_count();
669    let n_xyz = n_atoms * 3;
670    if coords.len() != n_xyz {
671        return Err(format!(
672            "coords length {} != 3 * atoms {}",
673            coords.len(),
674            n_atoms
675        ));
676    }
677    let elements: Vec<u8> = (0..n_atoms)
678        .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
679        .collect();
680    let mut grad = vec![0.0f64; n_xyz];
681    neb_energy_and_gradient(backend, smiles, &elements, &mol, coords, &mut grad)
682}
683
684/// Compute energy and return both energy (kcal/mol) and flat gradient (kcal/mol/Å).
685pub fn neb_backend_energy_and_gradient(
686    method: &str,
687    smiles: &str,
688    coords: &[f64],
689) -> Result<(f64, Vec<f64>), String> {
690    let backend = NebBackend::from_method(method)?;
691    let mol = crate::graph::Molecule::from_smiles(smiles)?;
692    let n_atoms = mol.graph.node_count();
693    let n_xyz = n_atoms * 3;
694    if coords.len() != n_xyz {
695        return Err(format!(
696            "coords length {} != 3 * atoms {}",
697            coords.len(),
698            n_atoms
699        ));
700    }
701    let elements: Vec<u8> = (0..n_atoms)
702        .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
703        .collect();
704    let mut grad = vec![0.0f64; n_xyz];
705    let energy = neb_energy_and_gradient(backend, smiles, &elements, &mol, coords, &mut grad)?;
706    Ok((energy, grad))
707}
708
709/// Compute energy and gradients using the specified backend.
710///
711/// Returns (energy_kcal_mol, gradients) or an error.
712pub fn compute_backend_energy_and_gradients(
713    backend: MdBackend,
714    smiles: &str,
715    elements: &[u8],
716    coords: &[f64],
717    grad: &mut [f64],
718) -> Result<f64, String> {
719    match backend {
720        MdBackend::Uff => {
721            let mol = crate::graph::Molecule::from_smiles(smiles)?;
722            let ff = crate::forcefield::builder::build_uff_force_field(&mol);
723            let e = ff.compute_system_energy_and_gradients(coords, grad);
724            Ok(e)
725        }
726        MdBackend::Pm3 => {
727            let positions: Vec<[f64; 3]> = coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
728            let grad_result = crate::pm3::gradients::compute_pm3_gradient(elements, &positions)?;
729            // PM3 gradient returns eV/Å; convert to kcal/mol/Å
730            let energy_kcal = grad_result.energy * 23.0605;
731            for (a, g) in grad_result.gradients.iter().enumerate() {
732                for d in 0..3 {
733                    grad[a * 3 + d] = g[d] * 23.0605;
734                }
735            }
736            Ok(energy_kcal)
737        }
738        MdBackend::Xtb => {
739            let positions: Vec<[f64; 3]> = coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
740            let grad_result = crate::xtb::gradients::compute_xtb_gradient(elements, &positions)?;
741            let energy_kcal = grad_result.energy * 23.0605;
742            for (a, g) in grad_result.gradients.iter().enumerate() {
743                for d in 0..3 {
744                    grad[a * 3 + d] = g[d] * 23.0605;
745                }
746            }
747            Ok(energy_kcal)
748        }
749    }
750}
751
752/// Run molecular dynamics using Velocity Verlet with a configurable force backend.
753///
754/// Supports UFF, PM3, and xTB backends. Note that PM3/xTB use numerical gradients
755/// and are significantly slower than UFF.
756pub fn simulate_velocity_verlet(
757    smiles: &str,
758    coords: &[f64],
759    elements: &[u8],
760    n_steps: usize,
761    dt_fs: f64,
762    seed: u64,
763    target_temp_and_tau: Option<(f64, f64)>,
764    backend: MdBackend,
765) -> Result<MdTrajectory, String> {
766    if n_steps == 0 {
767        return Err("n_steps must be > 0".to_string());
768    }
769    let n_atoms = elements.len();
770    if coords.len() != n_atoms * 3 {
771        return Err(format!(
772            "coords length {} != 3*atoms {}",
773            coords.len(),
774            n_atoms
775        ));
776    }
777
778    let masses_amu: Vec<f64> = elements.iter().map(|&z| atomic_mass_amu(z)).collect();
779
780    let mut x = coords.to_vec();
781    let mut grad = vec![0.0f64; n_atoms * 3];
782    let mut potential =
783        compute_backend_energy_and_gradients(backend, smiles, elements, &x, &mut grad)?;
784
785    let mut rng = StdRng::seed_from_u64(seed);
786    let mut v = vec![0.0f64; n_atoms * 3];
787
788    if let Some((target_temp_k, _)) = target_temp_and_tau {
789        for i in 0..n_atoms {
790            let sigma = ((R_GAS_KCAL_MOLK * target_temp_k)
791                / (masses_amu[i] * AMU_ANGFS2_TO_KCAL_MOL))
792                .sqrt();
793            v[3 * i] = sigma * sample_standard_normal(&mut rng);
794            v[3 * i + 1] = sigma * sample_standard_normal(&mut rng);
795            v[3 * i + 2] = sigma * sample_standard_normal(&mut rng);
796        }
797    }
798
799    let (ke0, t0) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
800    let mut frames = vec![MdFrame {
801        step: 0,
802        time_fs: 0.0,
803        coords: x.clone(),
804        potential_energy_kcal_mol: potential,
805        kinetic_energy_kcal_mol: ke0,
806        temperature_k: t0,
807    }];
808
809    for step in 1..=n_steps {
810        let mut a = vec![0.0f64; n_atoms * 3];
811        for i in 0..n_atoms {
812            let m = masses_amu[i];
813            for k in 0..3 {
814                a[3 * i + k] = -grad[3 * i + k] / (m * AMU_ANGFS2_TO_KCAL_MOL);
815            }
816        }
817
818        for i in 0..(n_atoms * 3) {
819            x[i] += v[i] * dt_fs + 0.5 * a[i] * dt_fs * dt_fs;
820        }
821
822        potential = compute_backend_energy_and_gradients(backend, smiles, elements, &x, &mut grad)?;
823
824        let mut a_new = vec![0.0f64; n_atoms * 3];
825        for i in 0..n_atoms {
826            let m = masses_amu[i];
827            for k in 0..3 {
828                a_new[3 * i + k] = -grad[3 * i + k] / (m * AMU_ANGFS2_TO_KCAL_MOL);
829            }
830        }
831
832        for i in 0..(n_atoms * 3) {
833            v[i] += 0.5 * (a[i] + a_new[i]) * dt_fs;
834        }
835
836        let (mut ke, mut temp_k) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
837        if let Some((target_temp_k, tau_fs)) = target_temp_and_tau {
838            let tau = tau_fs.max(1e-6);
839            let lambda = (1.0 + (dt_fs / tau) * (target_temp_k / temp_k.max(1e-6) - 1.0)).sqrt();
840            let lambda = lambda.clamp(0.5, 2.0);
841            for vi in &mut v {
842                *vi *= lambda;
843            }
844            let kt = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
845            ke = kt.0;
846            temp_k = kt.1;
847        }
848
849        if !x.iter().all(|v| v.is_finite()) || !potential.is_finite() {
850            return Err(format!("MD diverged at step {}", step));
851        }
852
853        frames.push(MdFrame {
854            step,
855            time_fs: step as f64 * dt_fs,
856            coords: x.clone(),
857            potential_energy_kcal_mol: potential,
858            kinetic_energy_kcal_mol: ke,
859            temperature_k: temp_k,
860        });
861    }
862
863    let energy_drift_percent = if target_temp_and_tau.is_none() {
864        let e0 = frames[0].potential_energy_kcal_mol + frames[0].kinetic_energy_kcal_mol;
865        let ef = frames
866            .last()
867            .map(|f| f.potential_energy_kcal_mol + f.kinetic_energy_kcal_mol)
868            .unwrap_or(e0);
869        if e0.abs() > 1e-10 {
870            Some((ef - e0) / e0.abs() * 100.0)
871        } else {
872            None
873        }
874    } else {
875        None
876    };
877
878    Ok(MdTrajectory {
879        frames,
880        dt_fs,
881        notes: vec![format!("Velocity Verlet with {:?} backend.", backend)],
882        energy_drift_percent,
883    })
884}
885
886/// Nosé-Hoover chain thermostat for rigorous NVT sampling.
887///
888/// Implements a chain of `chain_length` thermostats coupled to velocities.
889/// Produces canonical (NVT) ensemble with correct fluctuations.
890pub fn simulate_nose_hoover(
891    smiles: &str,
892    coords: &[f64],
893    elements: &[u8],
894    n_steps: usize,
895    dt_fs: f64,
896    target_temp_k: f64,
897    thermostat_mass: f64,
898    chain_length: usize,
899    seed: u64,
900    backend: MdBackend,
901) -> Result<MdTrajectory, String> {
902    if n_steps == 0 {
903        return Err("n_steps must be > 0".to_string());
904    }
905    let n_atoms = elements.len();
906    if coords.len() != n_atoms * 3 {
907        return Err("coords length mismatch".to_string());
908    }
909
910    let masses_amu: Vec<f64> = elements.iter().map(|&z| atomic_mass_amu(z)).collect();
911    let dof = (3 * n_atoms).saturating_sub(6).max(1) as f64;
912    let target_ke = 0.5 * dof * R_GAS_KCAL_MOLK * target_temp_k;
913
914    let mut x = coords.to_vec();
915    let mut grad = vec![0.0f64; n_atoms * 3];
916    let mut potential =
917        compute_backend_energy_and_gradients(backend, smiles, elements, &x, &mut grad)?;
918
919    let mut rng = StdRng::seed_from_u64(seed);
920    let mut v = vec![0.0f64; n_atoms * 3];
921    for i in 0..n_atoms {
922        let sigma =
923            ((R_GAS_KCAL_MOLK * target_temp_k) / (masses_amu[i] * AMU_ANGFS2_TO_KCAL_MOL)).sqrt();
924        v[3 * i] = sigma * sample_standard_normal(&mut rng);
925        v[3 * i + 1] = sigma * sample_standard_normal(&mut rng);
926        v[3 * i + 2] = sigma * sample_standard_normal(&mut rng);
927    }
928
929    // Nosé-Hoover chain variables
930    let chain_len = chain_length.max(1);
931    let q = vec![thermostat_mass; chain_len]; // thermostat masses
932    let mut xi = vec![0.0f64; chain_len]; // thermostat positions
933    let mut v_xi = vec![0.0f64; chain_len]; // thermostat velocities
934
935    let (ke0, t0) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
936    let mut frames = vec![MdFrame {
937        step: 0,
938        time_fs: 0.0,
939        coords: x.clone(),
940        potential_energy_kcal_mol: potential,
941        kinetic_energy_kcal_mol: ke0,
942        temperature_k: t0,
943    }];
944
945    let dt2 = dt_fs * 0.5;
946    let dt4 = dt_fs * 0.25;
947
948    for step in 1..=n_steps {
949        let (ke, _) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
950
951        // Nosé-Hoover chain: propagate thermostat (Yoshida-Suzuki-like)
952        // Force on first thermostat
953        let mut g_xi = vec![0.0f64; chain_len];
954        g_xi[0] = (2.0 * ke - 2.0 * target_ke) / q[0];
955        for j in 1..chain_len {
956            g_xi[j] =
957                (q[j - 1] * v_xi[j - 1] * v_xi[j - 1] - R_GAS_KCAL_MOLK * target_temp_k) / q[j];
958        }
959
960        // Update thermostat velocities (outside-in)
961        if chain_len > 1 {
962            v_xi[chain_len - 1] += g_xi[chain_len - 1] * dt4;
963        }
964        for j in (0..chain_len.saturating_sub(1)).rev() {
965            let exp_factor = (-v_xi[j + 1] * dt4).exp();
966            v_xi[j] = v_xi[j] * exp_factor + g_xi[j] * dt4;
967        }
968
969        // Scale velocities
970        let scale = (-v_xi[0] * dt2).exp();
971        for vi in v.iter_mut() {
972            *vi *= scale;
973        }
974
975        // Update thermostat positions
976        for j in 0..chain_len {
977            xi[j] += v_xi[j] * dt2;
978        }
979
980        // Recompute KE after scaling
981        let (ke_scaled, _) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
982        g_xi[0] = (2.0 * ke_scaled - 2.0 * target_ke) / q[0];
983
984        // Update thermostat velocities again (inside-out)
985        for j in 0..chain_len.saturating_sub(1) {
986            let exp_factor = (-v_xi[j + 1] * dt4).exp();
987            v_xi[j] = v_xi[j] * exp_factor + g_xi[j] * dt4;
988        }
989        if chain_len > 1 {
990            // Recompute force on last thermostat
991            g_xi[chain_len - 1] = (q[chain_len - 2] * v_xi[chain_len - 2] * v_xi[chain_len - 2]
992                - R_GAS_KCAL_MOLK * target_temp_k)
993                / q[chain_len - 1];
994            v_xi[chain_len - 1] += g_xi[chain_len - 1] * dt4;
995        }
996
997        // Velocity Verlet: half-step velocity
998        for i in 0..n_atoms {
999            let m = masses_amu[i];
1000            for k in 0..3 {
1001                v[3 * i + k] -= 0.5 * dt_fs * grad[3 * i + k] / (m * AMU_ANGFS2_TO_KCAL_MOL);
1002            }
1003        }
1004
1005        // Position update
1006        for i in 0..(n_atoms * 3) {
1007            x[i] += v[i] * dt_fs;
1008        }
1009
1010        // New forces
1011        potential = compute_backend_energy_and_gradients(backend, smiles, elements, &x, &mut grad)?;
1012
1013        // Velocity Verlet: second half-step velocity
1014        for i in 0..n_atoms {
1015            let m = masses_amu[i];
1016            for k in 0..3 {
1017                v[3 * i + k] -= 0.5 * dt_fs * grad[3 * i + k] / (m * AMU_ANGFS2_TO_KCAL_MOL);
1018            }
1019        }
1020
1021        let (ke_final, temp_k) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
1022
1023        if !x.iter().all(|v| v.is_finite()) || !potential.is_finite() {
1024            return Err(format!("NH-MD diverged at step {}", step));
1025        }
1026
1027        frames.push(MdFrame {
1028            step,
1029            time_fs: step as f64 * dt_fs,
1030            coords: x.clone(),
1031            potential_energy_kcal_mol: potential,
1032            kinetic_energy_kcal_mol: ke_final,
1033            temperature_k: temp_k,
1034        });
1035    }
1036
1037    // Temperature analysis
1038    let temps: Vec<f64> = frames
1039        .iter()
1040        .skip(frames.len() / 5)
1041        .map(|f| f.temperature_k)
1042        .collect();
1043    let avg_temp = temps.iter().sum::<f64>() / temps.len() as f64;
1044    let temp_variance =
1045        temps.iter().map(|t| (t - avg_temp).powi(2)).sum::<f64>() / temps.len() as f64;
1046
1047    Ok(MdTrajectory {
1048        frames,
1049        dt_fs,
1050        notes: vec![
1051            format!(
1052                "Nosé-Hoover chain thermostat (chain={}) with {:?} backend.",
1053                chain_len, backend
1054            ),
1055            format!(
1056                "Target T = {:.1} K, avg T = {:.1} K, σ(T) = {:.1} K",
1057                target_temp_k,
1058                avg_temp,
1059                temp_variance.sqrt()
1060            ),
1061        ],
1062        energy_drift_percent: None,
1063    })
1064}
1065
1066/// Export an MD trajectory to XYZ format string.
1067///
1068/// Each frame becomes a block with atom count, comment line (step/energy/temp),
1069/// and atom coordinates.
1070pub fn trajectory_to_xyz(trajectory: &MdTrajectory, elements: &[u8]) -> String {
1071    let n_atoms = elements.len();
1072    let mut output = String::new();
1073
1074    for frame in &trajectory.frames {
1075        output.push_str(&format!("{}\n", n_atoms));
1076        output.push_str(&format!(
1077            "step={} t={:.1}fs E_pot={:.4} E_kin={:.4} T={:.1}K\n",
1078            frame.step,
1079            frame.time_fs,
1080            frame.potential_energy_kcal_mol,
1081            frame.kinetic_energy_kcal_mol,
1082            frame.temperature_k,
1083        ));
1084        for i in 0..n_atoms {
1085            let sym = element_symbol(elements[i]);
1086            output.push_str(&format!(
1087                "{} {:.6} {:.6} {:.6}\n",
1088                sym,
1089                frame.coords[3 * i],
1090                frame.coords[3 * i + 1],
1091                frame.coords[3 * i + 2],
1092            ));
1093        }
1094    }
1095    output
1096}
1097
1098fn element_symbol(z: u8) -> &'static str {
1099    match z {
1100        1 => "H",
1101        2 => "He",
1102        3 => "Li",
1103        4 => "Be",
1104        5 => "B",
1105        6 => "C",
1106        7 => "N",
1107        8 => "O",
1108        9 => "F",
1109        10 => "Ne",
1110        11 => "Na",
1111        12 => "Mg",
1112        13 => "Al",
1113        14 => "Si",
1114        15 => "P",
1115        16 => "S",
1116        17 => "Cl",
1117        18 => "Ar",
1118        19 => "K",
1119        20 => "Ca",
1120        22 => "Ti",
1121        24 => "Cr",
1122        25 => "Mn",
1123        26 => "Fe",
1124        27 => "Co",
1125        28 => "Ni",
1126        29 => "Cu",
1127        30 => "Zn",
1128        35 => "Br",
1129        44 => "Ru",
1130        46 => "Pd",
1131        47 => "Ag",
1132        53 => "I",
1133        78 => "Pt",
1134        79 => "Au",
1135        _ => "X",
1136    }
1137}
1138
1139// ═══════════════════════════════════════════════════════════════════════════
1140//  Reaction dynamics: SMILES → embed → complex → NEB → full frame path
1141// ═══════════════════════════════════════════════════════════════════════════
1142
1143/// Configuration for a full reaction dynamics computation.
1144#[derive(Debug, Clone, Serialize, Deserialize)]
1145pub struct ReactionDynamicsConfig {
1146    /// Number of NEB images for the reactive region.
1147    pub n_neb_images: usize,
1148    /// NEB spring-coupled relaxation iterations.
1149    pub neb_iterations: usize,
1150    /// NEB spring constant (kcal/mol/Ų).
1151    pub spring_k: f64,
1152    /// NEB optimisation step size (Å).
1153    pub step_size: f64,
1154    /// Number of approach frames (molecules approaching).
1155    pub n_approach_frames: usize,
1156    /// Number of departure frames (products separating).
1157    pub n_departure_frames: usize,
1158    /// Far distance (Å) at start/end of approach/departure.
1159    pub far_distance: f64,
1160    /// Target distance between reactive atoms in the complexes (Å).
1161    pub reactive_distance: f64,
1162    /// Random seed for conformer generation.
1163    pub seed: u64,
1164}
1165
1166impl Default for ReactionDynamicsConfig {
1167    fn default() -> Self {
1168        Self {
1169            n_neb_images: 30,
1170            neb_iterations: 100,
1171            spring_k: 0.1,
1172            step_size: 0.01,
1173            n_approach_frames: 15,
1174            n_departure_frames: 15,
1175            far_distance: 8.0,
1176            reactive_distance: 2.0,
1177            seed: 42,
1178        }
1179    }
1180}
1181
1182/// A single frame along the reaction coordinate.
1183#[derive(Debug, Clone, Serialize, Deserialize)]
1184pub struct ReactionDynamicsFrame {
1185    /// Frame index (0-based).
1186    pub index: usize,
1187    /// Flat xyz coordinates `[x0,y0,z0, x1,y1,z1, ...]` in Å.
1188    pub coords: Vec<f64>,
1189    /// Energy at this frame (kcal/mol).
1190    pub energy_kcal_mol: f64,
1191    /// Phase label: `"approach"`, `"reaction"`, or `"departure"`.
1192    pub phase: String,
1193}
1194
1195/// Full result of a reaction dynamics computation.
1196#[derive(Debug, Clone, Serialize, Deserialize)]
1197pub struct ReactionDynamicsResult {
1198    /// All frames ordered along the reaction coordinate.
1199    pub frames: Vec<ReactionDynamicsFrame>,
1200    /// Atomic numbers for every atom in each frame (same for all frames).
1201    pub elements: Vec<u8>,
1202    /// Frame index of the transition state (highest energy).
1203    pub ts_frame_index: usize,
1204    /// Activation energy (kcal/mol) = E(TS) − E(first frame).
1205    pub activation_energy_kcal_mol: f64,
1206    /// Reaction energy (kcal/mol) = E(last frame) − E(first frame).
1207    pub reaction_energy_kcal_mol: f64,
1208    /// Method used for energy/gradient evaluation.
1209    pub method: String,
1210    /// Number of atoms per frame.
1211    pub n_atoms: usize,
1212    /// Informational notes.
1213    pub notes: Vec<String>,
1214}
1215
1216/// Centre-of-mass of a flat coordinate array.
1217pub fn com_flat(coords: &[f64]) -> [f64; 3] {
1218    let n = coords.len() / 3;
1219    if n == 0 {
1220        return [0.0; 3];
1221    }
1222    let mut cx = 0.0;
1223    let mut cy = 0.0;
1224    let mut cz = 0.0;
1225    for i in 0..n {
1226        cx += coords[i * 3];
1227        cy += coords[i * 3 + 1];
1228        cz += coords[i * 3 + 2];
1229    }
1230    let nf = n as f64;
1231    [cx / nf, cy / nf, cz / nf]
1232}
1233
1234/// Translate coords so COM is at origin.
1235pub fn centre_at_origin(coords: &mut [f64]) {
1236    let [cx, cy, cz] = com_flat(coords);
1237    for i in (0..coords.len()).step_by(3) {
1238        coords[i] -= cx;
1239        coords[i + 1] -= cy;
1240        coords[i + 2] -= cz;
1241    }
1242}
1243
1244/// Rodrigues rotation: rotate flat coords so direction `from` aligns with `to`.
1245fn rotate_to_align(coords: &mut [f64], from: [f64; 3], to: [f64; 3]) {
1246    let fl = (from[0] * from[0] + from[1] * from[1] + from[2] * from[2]).sqrt();
1247    let tl = (to[0] * to[0] + to[1] * to[1] + to[2] * to[2]).sqrt();
1248    if fl < 1e-10 || tl < 1e-10 {
1249        return;
1250    }
1251    let fx = from[0] / fl;
1252    let fy = from[1] / fl;
1253    let fz = from[2] / fl;
1254    let tx = to[0] / tl;
1255    let ty = to[1] / tl;
1256    let tz = to[2] / tl;
1257
1258    let kx = fy * tz - fz * ty;
1259    let ky = fz * tx - fx * tz;
1260    let kz = fx * ty - fy * tx;
1261    let sin_a = (kx * kx + ky * ky + kz * kz).sqrt();
1262    let cos_a = fx * tx + fy * ty + fz * tz;
1263
1264    if sin_a < 1e-10 {
1265        if cos_a > 0.0 {
1266            return; // already aligned
1267        }
1268        // Anti-aligned: reflect through a perpendicular plane
1269        for i in (0..coords.len()).step_by(3) {
1270            coords[i] = -coords[i];
1271        }
1272        return;
1273    }
1274
1275    let nkx = kx / sin_a;
1276    let nky = ky / sin_a;
1277    let nkz = kz / sin_a;
1278    let c = cos_a;
1279    let s = sin_a;
1280    let t1 = 1.0 - c;
1281
1282    for i in (0..coords.len()).step_by(3) {
1283        let x = coords[i];
1284        let y = coords[i + 1];
1285        let z = coords[i + 2];
1286        let dot = nkx * x + nky * y + nkz * z;
1287        coords[i] = x * c + (nky * z - nkz * y) * s + nkx * dot * t1;
1288        coords[i + 1] = y * c + (nkz * x - nkx * z) * s + nky * dot * t1;
1289        coords[i + 2] = z * c + (nkx * y - nky * x) * s + nkz * dot * t1;
1290    }
1291}
1292
1293/// Build a reactive complex from separately-embedded molecule conformers.
1294///
1295/// Each molecule is centred at its own COM, then the first two are oriented so
1296/// their closest-pair atoms face each other at `reactive_dist` Å apart.
1297/// Returns `(flat_coords, elements)` with the overall COM at origin.
1298fn build_reaction_complex(
1299    conformers: &[crate::ConformerResult],
1300    reactive_dist: f64,
1301) -> (Vec<f64>, Vec<u8>) {
1302    if conformers.is_empty() {
1303        return (vec![], vec![]);
1304    }
1305
1306    // Centre each molecule at its own COM
1307    let mut mols: Vec<(Vec<f64>, Vec<u8>)> = conformers
1308        .iter()
1309        .map(|c| {
1310            let mut coords = c.coords.clone();
1311            centre_at_origin(&mut coords);
1312            (coords, c.elements.clone())
1313        })
1314        .collect();
1315
1316    if mols.len() == 1 {
1317        let (mut coords, elems) = mols.remove(0);
1318        centre_at_origin(&mut coords);
1319        return (coords, elems);
1320    }
1321
1322    // Find closest inter-molecular pair between mol 0 and mol 1
1323    let n0 = mols[0].1.len();
1324    let n1 = mols[1].1.len();
1325    let mut best_i = 0usize;
1326    let mut best_j = 0usize;
1327    let mut best_d2 = f64::INFINITY;
1328    // Use a reference offset so we rank "face-to-face" proximity
1329    let ref_off = 4.0;
1330    for i in 0..n0 {
1331        let xi = mols[0].0[i * 3];
1332        let yi = mols[0].0[i * 3 + 1];
1333        let zi = mols[0].0[i * 3 + 2];
1334        for j in 0..n1 {
1335            let dx = mols[1].0[j * 3] + ref_off - xi;
1336            let dy = mols[1].0[j * 3 + 1] - yi;
1337            let dz = mols[1].0[j * 3 + 2] - zi;
1338            let d2 = dx * dx + dy * dy + dz * dz;
1339            if d2 < best_d2 {
1340                best_d2 = d2;
1341                best_i = i;
1342                best_j = j;
1343            }
1344        }
1345    }
1346
1347    // Orient mol 0 reactive atom → +X
1348    {
1349        let rx = mols[0].0[best_i * 3];
1350        let ry = mols[0].0[best_i * 3 + 1];
1351        let rz = mols[0].0[best_i * 3 + 2];
1352        let rl = (rx * rx + ry * ry + rz * rz).sqrt();
1353        if rl > 0.05 {
1354            rotate_to_align(&mut mols[0].0, [rx / rl, ry / rl, rz / rl], [1.0, 0.0, 0.0]);
1355        }
1356    }
1357
1358    // Orient mol 1 reactive atom → −X
1359    {
1360        let rx = mols[1].0[best_j * 3];
1361        let ry = mols[1].0[best_j * 3 + 1];
1362        let rz = mols[1].0[best_j * 3 + 2];
1363        let rl = (rx * rx + ry * ry + rz * rz).sqrt();
1364        if rl > 0.05 {
1365            rotate_to_align(
1366                &mut mols[1].0,
1367                [rx / rl, ry / rl, rz / rl],
1368                [-1.0, 0.0, 0.0],
1369            );
1370        }
1371    }
1372
1373    // Place mol 1 so reactive atoms are `reactive_dist` apart along X
1374    let ra_x = mols[0].0[best_i * 3];
1375    let rb_x = mols[1].0[best_j * 3];
1376    let offset_x = ra_x - rb_x + reactive_dist;
1377
1378    // Assemble combined coords
1379    let total_atoms: usize = mols.iter().map(|(_, e)| e.len()).sum();
1380    let mut all_coords = Vec::with_capacity(total_atoms * 3);
1381    let mut all_elements = Vec::with_capacity(total_atoms);
1382
1383    // Mol 0 — no offset
1384    all_coords.extend_from_slice(&mols[0].0);
1385    all_elements.extend_from_slice(&mols[0].1);
1386
1387    // Mol 1 — shifted along X
1388    for k in 0..n1 {
1389        all_coords.push(mols[1].0[k * 3] + offset_x);
1390        all_coords.push(mols[1].0[k * 3 + 1]);
1391        all_coords.push(mols[1].0[k * 3 + 2]);
1392    }
1393    all_elements.extend_from_slice(&mols[1].1);
1394
1395    // Spectator molecules (m ≥ 2) — stacked further along +X
1396    let mut extra = offset_x + 4.0;
1397    for mol in mols.iter().skip(2) {
1398        for k in 0..mol.1.len() {
1399            all_coords.push(mol.0[k * 3] + extra);
1400            all_coords.push(mol.0[k * 3 + 1]);
1401            all_coords.push(mol.0[k * 3 + 2]);
1402        }
1403        all_elements.extend_from_slice(&mol.1);
1404        extra += 4.0;
1405    }
1406
1407    centre_at_origin(&mut all_coords);
1408    (all_coords, all_elements)
1409}
1410
1411/// Greedy atom mapping between two complexes by element and distance.
1412///
1413/// Returns `mapping[i]` = index in the product complex that corresponds
1414/// to reactant atom `i`. Matches atoms of the same element, choosing the
1415/// closest unmatched pair greedily.
1416pub fn map_atoms_greedy(
1417    r_elements: &[u8],
1418    r_coords: &[f64],
1419    p_elements: &[u8],
1420    p_coords: &[f64],
1421) -> Result<Vec<usize>, String> {
1422    let n = r_elements.len();
1423    if n != p_elements.len() {
1424        return Err(format!(
1425            "Atom count mismatch: {} reactant atoms vs {} product atoms",
1426            n,
1427            p_elements.len(),
1428        ));
1429    }
1430
1431    // Group indices by element
1432    let mut r_by_elem: std::collections::HashMap<u8, Vec<usize>> = std::collections::HashMap::new();
1433    let mut p_by_elem: std::collections::HashMap<u8, Vec<usize>> = std::collections::HashMap::new();
1434    for i in 0..n {
1435        r_by_elem.entry(r_elements[i]).or_default().push(i);
1436        p_by_elem.entry(p_elements[i]).or_default().push(i);
1437    }
1438
1439    let mut mapping = vec![0usize; n];
1440    let mut used_p = vec![false; n];
1441
1442    for (&elem, r_indices) in &r_by_elem {
1443        let p_indices = p_by_elem
1444            .get(&elem)
1445            .ok_or_else(|| format!("Element Z={} present in reactants but not products", elem,))?;
1446        if r_indices.len() != p_indices.len() {
1447            return Err(format!(
1448                "Element Z={}: {} in reactants vs {} in products",
1449                elem,
1450                r_indices.len(),
1451                p_indices.len(),
1452            ));
1453        }
1454
1455        // Build distance pairs and sort by distance
1456        let mut pairs: Vec<(usize, usize, f64)> = Vec::new();
1457        for &ri in r_indices {
1458            for &pi in p_indices {
1459                let dx = r_coords[ri * 3] - p_coords[pi * 3];
1460                let dy = r_coords[ri * 3 + 1] - p_coords[pi * 3 + 1];
1461                let dz = r_coords[ri * 3 + 2] - p_coords[pi * 3 + 2];
1462                pairs.push((ri, pi, dx * dx + dy * dy + dz * dz));
1463            }
1464        }
1465        pairs.sort_by(|a, b| a.2.partial_cmp(&b.2).unwrap_or(std::cmp::Ordering::Equal));
1466
1467        let mut used_r = vec![false; n];
1468        for (ri, pi, _) in pairs {
1469            if used_r[ri] || used_p[pi] {
1470                continue;
1471            }
1472            mapping[ri] = pi;
1473            used_r[ri] = true;
1474            used_p[pi] = true;
1475        }
1476    }
1477
1478    Ok(mapping)
1479}
1480
1481/// Reorder product coords into reactant atom order using `mapping`.
1482pub fn reorder_coords(coords: &[f64], mapping: &[usize]) -> Vec<f64> {
1483    let mut out = vec![0.0; coords.len()];
1484    for (i, &j) in mapping.iter().enumerate() {
1485        out[i * 3] = coords[j * 3];
1486        out[i * 3 + 1] = coords[j * 3 + 1];
1487        out[i * 3 + 2] = coords[j * 3 + 2];
1488    }
1489    out
1490}
1491
1492/// Rigid-body slide: move each molecule's atoms toward/away from the global COM.
1493///
1494/// `mol_ranges[m] = (start_atom, end_atom)` (exclusive end) per molecule.
1495/// `alpha` = 0 → atoms at `far_dist` from global COM; `alpha` = 1 → at complex position.
1496pub fn slide_molecules(
1497    complex_coords: &[f64],
1498    mol_ranges: &[(usize, usize)],
1499    alpha: f64,
1500    far_dist: f64,
1501) -> Vec<f64> {
1502    let mut out = complex_coords.to_vec();
1503    let [gcx, gcy, gcz] = com_flat(complex_coords);
1504
1505    for &(start, end) in mol_ranges {
1506        let n = end - start;
1507        if n == 0 {
1508            continue;
1509        }
1510        let mut mx = 0.0;
1511        let mut my = 0.0;
1512        let mut mz = 0.0;
1513        for a in start..end {
1514            mx += complex_coords[a * 3];
1515            my += complex_coords[a * 3 + 1];
1516            mz += complex_coords[a * 3 + 2];
1517        }
1518        let nf = n as f64;
1519        mx /= nf;
1520        my /= nf;
1521        mz /= nf;
1522
1523        let dx = mx - gcx;
1524        let dy = my - gcy;
1525        let dz = mz - gcz;
1526        let d = (dx * dx + dy * dy + dz * dz).sqrt();
1527        if d < 0.01 {
1528            continue;
1529        }
1530        let target = d + (far_dist - d) * (1.0 - alpha);
1531        let shift = target - d;
1532        let nx = dx / d;
1533        let ny = dy / d;
1534        let nz = dz / d;
1535        for a in start..end {
1536            out[a * 3] += nx * shift;
1537            out[a * 3 + 1] += ny * shift;
1538            out[a * 3 + 2] += nz * shift;
1539        }
1540    }
1541    out
1542}
1543
1544/// Build a reactant complex guided by the product geometry.
1545///
1546/// For each reactant fragment, positions its embedded 3D geometry at the
1547/// centre-of-mass of the corresponding atoms in the **reordered product**
1548/// (product coords in reactant atom order).  Each fragment is Kabsch-aligned
1549/// to its product target positions and then translated to the product-derived
1550/// COM.  This ensures that the approach direction is chemically meaningful:
1551/// e.g. for H₂ + C₂H₂ → C₂H₄, H₂ approaches from beside the C≡C bond
1552/// (where the H atoms end up in ethylene), not end-on.
1553fn build_product_guided_reactant_complex(
1554    r_confs: &[crate::ConformerResult],
1555    p_reordered_coords: &[f64],
1556) -> Vec<f64> {
1557    let n_total: usize = r_confs.iter().map(|c| c.num_atoms).sum();
1558    let mut all_coords = vec![0.0f64; n_total * 3];
1559    let mut atom_off = 0usize;
1560
1561    for conf in r_confs {
1562        let n = conf.num_atoms;
1563        // Product-side positions for this fragment's atoms
1564        let p_frag: Vec<f64> = (atom_off..atom_off + n)
1565            .flat_map(|a| {
1566                [
1567                    p_reordered_coords[a * 3],
1568                    p_reordered_coords[a * 3 + 1],
1569                    p_reordered_coords[a * 3 + 2],
1570                ]
1571            })
1572            .collect();
1573        let p_com = com_flat(&p_frag);
1574
1575        // Reactant fragment's own geometry, centred at origin
1576        let mut r_frag = conf.coords.clone();
1577        centre_at_origin(&mut r_frag);
1578
1579        // Kabsch-align reactant fragment onto the product fragment positions
1580        // (both centred at their own COMs) to get the best rotation.
1581        if n >= 2 {
1582            let aligned = crate::alignment::kabsch::align_coordinates(&r_frag, &p_frag);
1583            // The aligned_coords are translated onto reference centroid,
1584            // so we must re-centre them and use the product COM.
1585            let ac = com_flat(&aligned.aligned_coords);
1586            for a in 0..n {
1587                all_coords[(atom_off + a) * 3] = aligned.aligned_coords[a * 3] - ac[0] + p_com[0];
1588                all_coords[(atom_off + a) * 3 + 1] =
1589                    aligned.aligned_coords[a * 3 + 1] - ac[1] + p_com[1];
1590                all_coords[(atom_off + a) * 3 + 2] =
1591                    aligned.aligned_coords[a * 3 + 2] - ac[2] + p_com[2];
1592            }
1593        } else {
1594            // Single atom — just place at product position
1595            for a in 0..n {
1596                all_coords[(atom_off + a) * 3] = r_frag[a * 3] + p_com[0];
1597                all_coords[(atom_off + a) * 3 + 1] = r_frag[a * 3 + 1] + p_com[1];
1598                all_coords[(atom_off + a) * 3 + 2] = r_frag[a * 3 + 2] + p_com[2];
1599            }
1600        }
1601
1602        atom_off += n;
1603    }
1604
1605    centre_at_origin(&mut all_coords);
1606    all_coords
1607}
1608
1609/// Compute a full reaction dynamics path: embed reactants + products,
1610/// build oriented complexes, run NEB for the reactive region, and generate
1611/// approach/departure frames — all energies computed in Rust with the
1612/// chosen quantum-chemistry method.
1613///
1614/// The reactant complex orientation is **derived from the product geometry**:
1615/// each reactant fragment is positioned where its atoms end up in the product,
1616/// ensuring a chemically meaningful approach direction.  For example, in
1617/// H₂ + C₂H₂ → C₂H₄ the H₂ approaches side-on to the π system rather than
1618/// end-on along the C≡C axis.
1619///
1620/// # Arguments
1621///
1622/// * `reactant_smiles` — SMILES of each reactant fragment (e.g. `["CC(=O)O", "N"]`).
1623/// * `product_smiles`  — SMILES of each product fragment (e.g. `["CC(=O)N", "O"]`).
1624/// * `method`          — NEB backend: `"uff"`, `"mmff94"`, `"pm3"`, `"xtb"`, `"gfn1"`, `"gfn2"`, `"hf3c"`.
1625/// * `config`          — [`ReactionDynamicsConfig`] with NEB/path parameters.
1626///
1627/// Returns [`ReactionDynamicsResult`] with all frames, energies, and TS info.
1628pub fn compute_reaction_dynamics(
1629    reactant_smiles: &[&str],
1630    product_smiles: &[&str],
1631    method: &str,
1632    config: &ReactionDynamicsConfig,
1633) -> Result<ReactionDynamicsResult, String> {
1634    if reactant_smiles.is_empty() {
1635        return Err("At least one reactant SMILES is required".into());
1636    }
1637    if product_smiles.is_empty() {
1638        return Err("At least one product SMILES is required".into());
1639    }
1640
1641    let backend = NebBackend::from_method(method)?;
1642
1643    // ── 1. Embed all fragments ─────────────────────────────────────────
1644    let r_confs: Vec<crate::ConformerResult> = reactant_smiles
1645        .iter()
1646        .map(|s| crate::embed(s, config.seed))
1647        .collect();
1648    for (i, c) in r_confs.iter().enumerate() {
1649        if let Some(ref e) = c.error {
1650            return Err(format!(
1651                "Failed to embed reactant '{}': {}",
1652                reactant_smiles[i], e
1653            ));
1654        }
1655    }
1656
1657    let p_confs: Vec<crate::ConformerResult> = product_smiles
1658        .iter()
1659        .map(|s| crate::embed(s, config.seed))
1660        .collect();
1661    for (i, c) in p_confs.iter().enumerate() {
1662        if let Some(ref e) = c.error {
1663            return Err(format!(
1664                "Failed to embed product '{}': {}",
1665                product_smiles[i], e
1666            ));
1667        }
1668    }
1669
1670    // Validate atom conservation
1671    let r_total: usize = r_confs.iter().map(|c| c.num_atoms).sum();
1672    let p_total: usize = p_confs.iter().map(|c| c.num_atoms).sum();
1673    if r_total != p_total {
1674        return Err(format!(
1675            "Atom count mismatch: {} atoms in reactants vs {} in products — \
1676             atoms must be conserved",
1677            r_total, p_total,
1678        ));
1679    }
1680
1681    // ── 2. Build product complex and collect elements ──────────────────
1682    let (p_coords, p_elements) = build_reaction_complex(&p_confs, config.reactive_distance);
1683
1684    // Reactant elements in concatenation order
1685    let r_elements: Vec<u8> = r_confs
1686        .iter()
1687        .flat_map(|c| c.elements.iter().copied())
1688        .collect();
1689
1690    // ── 3. Atom mapping: greedy by element + distance ──────────────────
1691    //  We need a temporary reactant complex just for the distance mapping.
1692    //  Use product-agnostic placement first, then rebuild with guidance.
1693    let (r_coords_tmp, _) = build_reaction_complex(&r_confs, config.reactive_distance);
1694    let mapping = map_atoms_greedy(&r_elements, &r_coords_tmp, &p_elements, &p_coords)?;
1695    let p_reordered = reorder_coords(&p_coords, &mapping);
1696
1697    // ── 4. Build reactant complex guided by product geometry ───────────
1698    let r_coords = build_product_guided_reactant_complex(&r_confs, &p_reordered);
1699
1700    // ── 5. Combined SMILES (dot-separated reactants) for NEB topology ──
1701    let combined_smiles = reactant_smiles.join(".");
1702
1703    // ── 6. NEB between product-guided reactant complex and product ─────
1704    let neb = compute_simplified_neb_path_configurable(
1705        &combined_smiles,
1706        &r_coords,
1707        &p_reordered,
1708        config.n_neb_images,
1709        config.neb_iterations,
1710        config.spring_k,
1711        config.step_size,
1712        method,
1713    )?;
1714
1715    // ── 7. Compute molecule ranges for approach sliding ────────────────
1716    let r_mol_ranges: Vec<(usize, usize)> = {
1717        let mut ranges = Vec::new();
1718        let mut off = 0usize;
1719        for c in &r_confs {
1720            ranges.push((off, off + c.num_atoms));
1721            off += c.num_atoms;
1722        }
1723        ranges
1724    };
1725
1726    // ── 8. Build approach frames with energies ─────────────────────────
1727    let mut frames = Vec::new();
1728    let mol = crate::graph::Molecule::from_smiles(&combined_smiles)?;
1729    let n_xyz = r_total * 3;
1730
1731    let na = config.n_approach_frames;
1732    for i in 0..na {
1733        let alpha = if na > 1 {
1734            i as f64 / (na - 1) as f64
1735        } else {
1736            1.0
1737        };
1738        let coords = slide_molecules(&r_coords, &r_mol_ranges, alpha, config.far_distance);
1739        let mut grad = vec![0.0; n_xyz];
1740        let energy = neb_energy_and_gradient(
1741            backend,
1742            &combined_smiles,
1743            &r_elements,
1744            &mol,
1745            &coords,
1746            &mut grad,
1747        )
1748        .unwrap_or(0.0);
1749        frames.push(ReactionDynamicsFrame {
1750            index: frames.len(),
1751            coords,
1752            energy_kcal_mol: energy,
1753            phase: "approach".into(),
1754        });
1755    }
1756
1757    // ── 9. Reaction frames from NEB ────────────────────────────────────
1758    for img in &neb.images {
1759        frames.push(ReactionDynamicsFrame {
1760            index: frames.len(),
1761            coords: img.coords.clone(),
1762            energy_kcal_mol: img.potential_energy_kcal_mol,
1763            phase: "reaction".into(),
1764        });
1765    }
1766
1767    // ── 10. Departure frames with energies ─────────────────────────────
1768    //   If all products form a single molecule (e.g. C2H2 + H2 → C2H4),
1769    //   the departure phase holds the product geometry steady — no sliding.
1770    //   When there are multiple product molecules, slide them apart.
1771    let nd = config.n_departure_frames;
1772    let single_product = p_confs.len() == 1;
1773
1774    // Product molecule ranges (after reordering into reactant atom order)
1775    let p_mol_ranges: Vec<(usize, usize)> = if single_product {
1776        vec![(0, r_total)]
1777    } else {
1778        let p_mol_assign: Vec<usize> = {
1779            let mut assign = vec![0usize; p_total];
1780            let mut off = 0usize;
1781            for (m, c) in p_confs.iter().enumerate() {
1782                for a in off..off + c.num_atoms {
1783                    assign[a] = m;
1784                }
1785                off += c.num_atoms;
1786            }
1787            mapping.iter().map(|&pi| assign[pi]).collect()
1788        };
1789        let n_mols = p_confs.len();
1790        let mut ranges = Vec::with_capacity(n_mols);
1791        for m in 0..n_mols {
1792            let atoms: Vec<usize> = p_mol_assign
1793                .iter()
1794                .enumerate()
1795                .filter(|(_, &a)| a == m)
1796                .map(|(i, _)| i)
1797                .collect();
1798            if let (Some(&lo), Some(&hi)) = (atoms.first(), atoms.last()) {
1799                ranges.push((lo, hi + 1));
1800            }
1801        }
1802        ranges
1803    };
1804
1805    for i in 0..nd {
1806        let alpha = if nd > 1 {
1807            1.0 - i as f64 / (nd - 1) as f64
1808        } else {
1809            0.0
1810        };
1811        let coords = if single_product {
1812            // Single product — no separation; hold the product geometry
1813            p_reordered.clone()
1814        } else {
1815            slide_molecules(&p_reordered, &p_mol_ranges, alpha, config.far_distance)
1816        };
1817        let mut grad = vec![0.0; n_xyz];
1818        let energy = neb_energy_and_gradient(
1819            backend,
1820            &combined_smiles,
1821            &r_elements,
1822            &mol,
1823            &coords,
1824            &mut grad,
1825        )
1826        .unwrap_or(0.0);
1827        frames.push(ReactionDynamicsFrame {
1828            index: frames.len(),
1829            coords,
1830            energy_kcal_mol: energy,
1831            phase: "departure".into(),
1832        });
1833    }
1834
1835    // ── 11. Find TS and compute energetics ─────────────────────────────
1836    let ts_idx = frames
1837        .iter()
1838        .enumerate()
1839        .max_by(|(_, a), (_, b)| {
1840            a.energy_kcal_mol
1841                .partial_cmp(&b.energy_kcal_mol)
1842                .unwrap_or(std::cmp::Ordering::Equal)
1843        })
1844        .map(|(i, _)| i)
1845        .unwrap_or(0);
1846
1847    let e_first = frames.first().map(|f| f.energy_kcal_mol).unwrap_or(0.0);
1848    let e_last = frames.last().map(|f| f.energy_kcal_mol).unwrap_or(0.0);
1849    let e_ts = frames[ts_idx].energy_kcal_mol;
1850
1851    let mut notes = neb.notes;
1852    notes.push(format!(
1853        "Full reaction path: {} approach + {} NEB + {} departure = {} total frames.",
1854        na,
1855        neb.images.len(),
1856        nd,
1857        frames.len(),
1858    ));
1859    notes.push(
1860        "Reactant complex oriented using product geometry for chemically meaningful approach."
1861            .to_string(),
1862    );
1863
1864    Ok(ReactionDynamicsResult {
1865        frames,
1866        elements: r_elements,
1867        ts_frame_index: ts_idx,
1868        activation_energy_kcal_mol: e_ts - e_first,
1869        reaction_energy_kcal_mol: e_last - e_first,
1870        method: method.to_string(),
1871        n_atoms: r_total,
1872        notes,
1873    })
1874}