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;
7
8/// Force backend for molecular dynamics.
9#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
10pub enum MdBackend {
11    /// UFF force field (fast, approximate).
12    Uff,
13    /// PM3 semi-empirical (slower, more accurate).
14    Pm3,
15    /// GFN-xTB tight-binding (moderate speed, good for metals).
16    Xtb,
17}
18
19/// One trajectory frame for molecular-dynamics sampling.
20#[derive(Debug, Clone, Serialize, Deserialize)]
21pub struct MdFrame {
22    /// Integration step index.
23    pub step: usize,
24    /// Elapsed simulation time in femtoseconds.
25    pub time_fs: f64,
26    /// Flat xyz coordinates in angstroms.
27    pub coords: Vec<f64>,
28    /// Potential energy from UFF in kcal/mol.
29    pub potential_energy_kcal_mol: f64,
30    /// Kinetic energy in kcal/mol.
31    pub kinetic_energy_kcal_mol: f64,
32    /// Instantaneous temperature estimate in K.
33    pub temperature_k: f64,
34}
35
36/// Full trajectory output for exploratory molecular dynamics.
37#[derive(Debug, Clone, Serialize, Deserialize)]
38pub struct MdTrajectory {
39    /// Stored MD frames.
40    pub frames: Vec<MdFrame>,
41    /// Timestep in femtoseconds.
42    pub dt_fs: f64,
43    /// Notes and caveats for interpretation.
44    pub notes: Vec<String>,
45    /// Energy conservation drift: (E_total_final - E_total_initial) / |E_total_initial| * 100%.
46    /// Only meaningful for NVE (no thermostat) simulations.
47    pub energy_drift_percent: Option<f64>,
48}
49
50/// One image (node) on a simplified NEB path.
51#[derive(Debug, Clone, Serialize, Deserialize)]
52pub struct NebImage {
53    /// Image index from reactant (0) to product (n-1).
54    pub index: usize,
55    /// Flat xyz coordinates in angstroms.
56    pub coords: Vec<f64>,
57    /// UFF potential energy in kcal/mol.
58    pub potential_energy_kcal_mol: f64,
59}
60
61/// Simplified NEB output for low-cost pathway exploration.
62#[derive(Debug, Clone, Serialize, Deserialize)]
63pub struct NebPathResult {
64    /// Ordered path images from reactant to product.
65    pub images: Vec<NebImage>,
66    /// Notes and caveats.
67    pub notes: Vec<String>,
68}
69
70pub fn atomic_mass_amu(z: u8) -> f64 {
71    match z {
72        1 => 1.008,
73        5 => 10.81,
74        6 => 12.011,
75        7 => 14.007,
76        8 => 15.999,
77        9 => 18.998,
78        14 => 28.085,
79        15 => 30.974,
80        16 => 32.06,
81        17 => 35.45,
82        35 => 79.904,
83        53 => 126.904,
84        26 => 55.845,
85        46 => 106.42,
86        78 => 195.084,
87        _ => 12.0,
88    }
89}
90
91fn sample_standard_normal(rng: &mut StdRng) -> f64 {
92    let u1 = (1.0 - rng.gen::<f64>()).max(1e-12);
93    let u2 = rng.gen::<f64>();
94    (-2.0 * u1.ln()).sqrt() * (2.0 * std::f64::consts::PI * u2).cos()
95}
96
97fn kinetic_energy_and_temperature(
98    velocities: &[f64],
99    masses_amu: &[f64],
100    n_atoms: usize,
101) -> (f64, f64) {
102    let mut ke = 0.0;
103    for i in 0..n_atoms {
104        let vx = velocities[3 * i];
105        let vy = velocities[3 * i + 1];
106        let vz = velocities[3 * i + 2];
107        let v2 = vx * vx + vy * vy + vz * vz;
108        ke += 0.5 * masses_amu[i] * v2 * AMU_ANGFS2_TO_KCAL_MOL;
109    }
110    let dof = (3 * n_atoms).saturating_sub(6).max(1) as f64;
111    let t = 2.0 * ke / (dof * R_GAS_KCAL_MOLK);
112    (ke, t)
113}
114
115/// Run short exploratory molecular dynamics using Velocity Verlet and optional Berendsen NVT.
116pub fn simulate_velocity_verlet_uff(
117    smiles: &str,
118    coords: &[f64],
119    n_steps: usize,
120    dt_fs: f64,
121    seed: u64,
122    target_temp_and_tau: Option<(f64, f64)>,
123) -> Result<MdTrajectory, String> {
124    if n_steps == 0 {
125        return Err("n_steps must be > 0".to_string());
126    }
127    if dt_fs <= 0.0 {
128        return Err("dt_fs must be > 0".to_string());
129    }
130
131    let mol = crate::graph::Molecule::from_smiles(smiles)?;
132    let n_atoms = mol.graph.node_count();
133    if coords.len() != n_atoms * 3 {
134        return Err(format!(
135            "coords length {} != 3 * atoms {}",
136            coords.len(),
137            n_atoms
138        ));
139    }
140
141    let masses_amu: Vec<f64> = (0..n_atoms)
142        .map(petgraph::graph::NodeIndex::new)
143        .map(|idx| atomic_mass_amu(mol.graph[idx].element))
144        .collect();
145
146    let ff = crate::forcefield::builder::build_uff_force_field(&mol);
147    let mut x = coords.to_vec();
148    let mut grad = vec![0.0f64; n_atoms * 3];
149    let mut potential = ff.compute_system_energy_and_gradients(&x, &mut grad);
150
151    let mut rng = StdRng::seed_from_u64(seed);
152    let mut v = vec![0.0f64; n_atoms * 3];
153
154    if let Some((target_temp_k, _tau_fs)) = target_temp_and_tau {
155        for i in 0..n_atoms {
156            let sigma = ((R_GAS_KCAL_MOLK * target_temp_k)
157                / (masses_amu[i] * AMU_ANGFS2_TO_KCAL_MOL))
158                .sqrt();
159            v[3 * i] = sigma * sample_standard_normal(&mut rng);
160            v[3 * i + 1] = sigma * sample_standard_normal(&mut rng);
161            v[3 * i + 2] = sigma * sample_standard_normal(&mut rng);
162        }
163    }
164
165    let (ke0, t0) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
166    let mut frames = vec![MdFrame {
167        step: 0,
168        time_fs: 0.0,
169        coords: x.clone(),
170        potential_energy_kcal_mol: potential,
171        kinetic_energy_kcal_mol: ke0,
172        temperature_k: t0,
173    }];
174
175    for step in 1..=n_steps {
176        let mut a = vec![0.0f64; n_atoms * 3];
177        for i in 0..n_atoms {
178            let m = masses_amu[i];
179            a[3 * i] = -grad[3 * i] / (m * AMU_ANGFS2_TO_KCAL_MOL);
180            a[3 * i + 1] = -grad[3 * i + 1] / (m * AMU_ANGFS2_TO_KCAL_MOL);
181            a[3 * i + 2] = -grad[3 * i + 2] / (m * AMU_ANGFS2_TO_KCAL_MOL);
182        }
183
184        for i in 0..(n_atoms * 3) {
185            x[i] += v[i] * dt_fs + 0.5 * a[i] * dt_fs * dt_fs;
186        }
187
188        potential = ff.compute_system_energy_and_gradients(&x, &mut grad);
189
190        let mut a_new = vec![0.0f64; n_atoms * 3];
191        for i in 0..n_atoms {
192            let m = masses_amu[i];
193            a_new[3 * i] = -grad[3 * i] / (m * AMU_ANGFS2_TO_KCAL_MOL);
194            a_new[3 * i + 1] = -grad[3 * i + 1] / (m * AMU_ANGFS2_TO_KCAL_MOL);
195            a_new[3 * i + 2] = -grad[3 * i + 2] / (m * AMU_ANGFS2_TO_KCAL_MOL);
196        }
197
198        for i in 0..(n_atoms * 3) {
199            v[i] += 0.5 * (a[i] + a_new[i]) * dt_fs;
200        }
201
202        let (mut ke, mut temp_k) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
203        if let Some((target_temp_k, tau_fs)) = target_temp_and_tau {
204            let tau = tau_fs.max(1e-6);
205            let lambda = (1.0 + (dt_fs / tau) * (target_temp_k / temp_k.max(1e-6) - 1.0)).sqrt();
206            let lambda = lambda.clamp(0.5, 2.0);
207            for vi in &mut v {
208                *vi *= lambda;
209            }
210            let kt = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
211            ke = kt.0;
212            temp_k = kt.1;
213        }
214
215        if !x.iter().all(|v| v.is_finite()) || !potential.is_finite() || !ke.is_finite() {
216            return Err(format!(
217                "MD diverged at step {} (non-finite coordinates/energy)",
218                step
219            ));
220        }
221
222        frames.push(MdFrame {
223            step,
224            time_fs: step as f64 * dt_fs,
225            coords: x.clone(),
226            potential_energy_kcal_mol: potential,
227            kinetic_energy_kcal_mol: ke,
228            temperature_k: temp_k,
229        });
230    }
231
232    let mut notes = vec![
233        "Velocity Verlet integration over UFF force-field gradients for short exploratory trajectories."
234            .to_string(),
235    ];
236    if target_temp_and_tau.is_some() {
237        notes.push(
238            "Berendsen thermostat rescaling applied for approximate constant-temperature sampling."
239                .to_string(),
240        );
241    } else {
242        notes.push(
243            "No thermostat applied (NVE-like propagation with current numerical approximations)."
244                .to_string(),
245        );
246    }
247
248    Ok(MdTrajectory {
249        energy_drift_percent: if target_temp_and_tau.is_none() {
250            let e0 = frames[0].potential_energy_kcal_mol + frames[0].kinetic_energy_kcal_mol;
251            let ef = frames
252                .last()
253                .map(|f| f.potential_energy_kcal_mol + f.kinetic_energy_kcal_mol)
254                .unwrap_or(e0);
255            if e0.abs() > 1e-10 {
256                Some((ef - e0) / e0.abs() * 100.0)
257            } else {
258                None
259            }
260        } else {
261            None
262        },
263        frames,
264        dt_fs,
265        notes,
266    })
267}
268
269/// Build a simplified NEB-like path using linear interpolation and spring-coupled relaxation.
270pub fn compute_simplified_neb_path(
271    smiles: &str,
272    start_coords: &[f64],
273    end_coords: &[f64],
274    n_images: usize,
275    n_iter: usize,
276    spring_k: f64,
277    step_size: f64,
278) -> Result<NebPathResult, String> {
279    if n_images < 2 {
280        return Err("n_images must be >= 2".to_string());
281    }
282    if n_iter == 0 {
283        return Err("n_iter must be > 0".to_string());
284    }
285    if step_size <= 0.0 {
286        return Err("step_size must be > 0".to_string());
287    }
288
289    let mol = crate::graph::Molecule::from_smiles(smiles)?;
290    let n_atoms = mol.graph.node_count();
291    let n_xyz = n_atoms * 3;
292    if start_coords.len() != n_xyz || end_coords.len() != n_xyz {
293        return Err(format!(
294            "start/end coords must each have length {} (3 * n_atoms)",
295            n_xyz
296        ));
297    }
298
299    let ff = crate::forcefield::builder::build_uff_force_field(&mol);
300
301    let mut images = vec![vec![0.0f64; n_xyz]; n_images];
302    for (img_idx, img) in images.iter_mut().enumerate() {
303        let t = img_idx as f64 / (n_images - 1) as f64;
304        for k in 0..n_xyz {
305            img[k] = (1.0 - t) * start_coords[k] + t * end_coords[k];
306        }
307    }
308
309    for _ in 0..n_iter {
310        let prev = images.clone();
311
312        #[cfg(feature = "parallel")]
313        {
314            use rayon::prelude::*;
315            let updated: Vec<(usize, Vec<f64>)> = (1..(n_images - 1))
316                .into_par_iter()
317                .map(|i| {
318                    let ff_local = crate::forcefield::builder::build_uff_force_field(&mol);
319                    let mut grad = vec![0.0f64; n_xyz];
320                    let _ = ff_local.compute_system_energy_and_gradients(&prev[i], &mut grad);
321                    let mut new_img = prev[i].clone();
322                    for k in 0..n_xyz {
323                        let spring_force =
324                            spring_k * (prev[i + 1][k] - 2.0 * prev[i][k] + prev[i - 1][k]);
325                        let total_force = -grad[k] + spring_force;
326                        new_img[k] = prev[i][k] + step_size * total_force;
327                    }
328                    (i, new_img)
329                })
330                .collect();
331            for (i, img) in updated {
332                images[i] = img;
333            }
334        }
335
336        #[cfg(not(feature = "parallel"))]
337        {
338            for i in 1..(n_images - 1) {
339                let mut grad = vec![0.0f64; n_xyz];
340                let _ = ff.compute_system_energy_and_gradients(&prev[i], &mut grad);
341                for k in 0..n_xyz {
342                    let spring_force =
343                        spring_k * (prev[i + 1][k] - 2.0 * prev[i][k] + prev[i - 1][k]);
344                    let total_force = -grad[k] + spring_force;
345                    images[i][k] = prev[i][k] + step_size * total_force;
346                }
347            }
348        }
349    }
350
351    let mut out_images = Vec::with_capacity(n_images);
352    for (i, coords) in images.into_iter().enumerate() {
353        let mut grad = vec![0.0f64; n_xyz];
354        let e = ff.compute_system_energy_and_gradients(&coords, &mut grad);
355        out_images.push(NebImage {
356            index: i,
357            coords,
358            potential_energy_kcal_mol: e,
359        });
360    }
361
362    Ok(NebPathResult {
363        images: out_images,
364        notes: vec![
365            "Simplified NEB: linear interpolation + spring-coupled UFF gradient relaxation on internal images."
366                .to_string(),
367            "This is a low-cost exploratory path tool and not a full climbing-image / tangent-projected NEB implementation."
368                .to_string(),
369        ],
370    })
371}
372
373/// Compute energy and gradients using the specified backend.
374///
375/// Returns (energy_kcal_mol, gradients) or an error.
376pub fn compute_backend_energy_and_gradients(
377    backend: MdBackend,
378    smiles: &str,
379    elements: &[u8],
380    coords: &[f64],
381    grad: &mut [f64],
382) -> Result<f64, String> {
383    match backend {
384        MdBackend::Uff => {
385            let mol = crate::graph::Molecule::from_smiles(smiles)?;
386            let ff = crate::forcefield::builder::build_uff_force_field(&mol);
387            let e = ff.compute_system_energy_and_gradients(coords, grad);
388            Ok(e)
389        }
390        MdBackend::Pm3 => {
391            let positions: Vec<[f64; 3]> = coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
392            let grad_result = crate::pm3::gradients::compute_pm3_gradient(elements, &positions)?;
393            // PM3 gradient returns eV/Å; convert to kcal/mol/Å
394            let energy_kcal = grad_result.energy * 23.0605;
395            for (a, g) in grad_result.gradients.iter().enumerate() {
396                for d in 0..3 {
397                    grad[a * 3 + d] = g[d] * 23.0605;
398                }
399            }
400            Ok(energy_kcal)
401        }
402        MdBackend::Xtb => {
403            let positions: Vec<[f64; 3]> = coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
404            let grad_result = crate::xtb::gradients::compute_xtb_gradient(elements, &positions)?;
405            let energy_kcal = grad_result.energy * 23.0605;
406            for (a, g) in grad_result.gradients.iter().enumerate() {
407                for d in 0..3 {
408                    grad[a * 3 + d] = g[d] * 23.0605;
409                }
410            }
411            Ok(energy_kcal)
412        }
413    }
414}
415
416/// Run molecular dynamics using Velocity Verlet with a configurable force backend.
417///
418/// Supports UFF, PM3, and xTB backends. Note that PM3/xTB use numerical gradients
419/// and are significantly slower than UFF.
420pub fn simulate_velocity_verlet(
421    smiles: &str,
422    coords: &[f64],
423    elements: &[u8],
424    n_steps: usize,
425    dt_fs: f64,
426    seed: u64,
427    target_temp_and_tau: Option<(f64, f64)>,
428    backend: MdBackend,
429) -> Result<MdTrajectory, String> {
430    if n_steps == 0 {
431        return Err("n_steps must be > 0".to_string());
432    }
433    let n_atoms = elements.len();
434    if coords.len() != n_atoms * 3 {
435        return Err(format!(
436            "coords length {} != 3*atoms {}",
437            coords.len(),
438            n_atoms
439        ));
440    }
441
442    let masses_amu: Vec<f64> = elements.iter().map(|&z| atomic_mass_amu(z)).collect();
443
444    let mut x = coords.to_vec();
445    let mut grad = vec![0.0f64; n_atoms * 3];
446    let mut potential =
447        compute_backend_energy_and_gradients(backend, smiles, elements, &x, &mut grad)?;
448
449    let mut rng = StdRng::seed_from_u64(seed);
450    let mut v = vec![0.0f64; n_atoms * 3];
451
452    if let Some((target_temp_k, _)) = target_temp_and_tau {
453        for i in 0..n_atoms {
454            let sigma = ((R_GAS_KCAL_MOLK * target_temp_k)
455                / (masses_amu[i] * AMU_ANGFS2_TO_KCAL_MOL))
456                .sqrt();
457            v[3 * i] = sigma * sample_standard_normal(&mut rng);
458            v[3 * i + 1] = sigma * sample_standard_normal(&mut rng);
459            v[3 * i + 2] = sigma * sample_standard_normal(&mut rng);
460        }
461    }
462
463    let (ke0, t0) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
464    let mut frames = vec![MdFrame {
465        step: 0,
466        time_fs: 0.0,
467        coords: x.clone(),
468        potential_energy_kcal_mol: potential,
469        kinetic_energy_kcal_mol: ke0,
470        temperature_k: t0,
471    }];
472
473    for step in 1..=n_steps {
474        let mut a = vec![0.0f64; n_atoms * 3];
475        for i in 0..n_atoms {
476            let m = masses_amu[i];
477            for k in 0..3 {
478                a[3 * i + k] = -grad[3 * i + k] / (m * AMU_ANGFS2_TO_KCAL_MOL);
479            }
480        }
481
482        for i in 0..(n_atoms * 3) {
483            x[i] += v[i] * dt_fs + 0.5 * a[i] * dt_fs * dt_fs;
484        }
485
486        potential = compute_backend_energy_and_gradients(backend, smiles, elements, &x, &mut grad)?;
487
488        let mut a_new = vec![0.0f64; n_atoms * 3];
489        for i in 0..n_atoms {
490            let m = masses_amu[i];
491            for k in 0..3 {
492                a_new[3 * i + k] = -grad[3 * i + k] / (m * AMU_ANGFS2_TO_KCAL_MOL);
493            }
494        }
495
496        for i in 0..(n_atoms * 3) {
497            v[i] += 0.5 * (a[i] + a_new[i]) * dt_fs;
498        }
499
500        let (mut ke, mut temp_k) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
501        if let Some((target_temp_k, tau_fs)) = target_temp_and_tau {
502            let tau = tau_fs.max(1e-6);
503            let lambda = (1.0 + (dt_fs / tau) * (target_temp_k / temp_k.max(1e-6) - 1.0)).sqrt();
504            let lambda = lambda.clamp(0.5, 2.0);
505            for vi in &mut v {
506                *vi *= lambda;
507            }
508            let kt = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
509            ke = kt.0;
510            temp_k = kt.1;
511        }
512
513        if !x.iter().all(|v| v.is_finite()) || !potential.is_finite() {
514            return Err(format!("MD diverged at step {}", step));
515        }
516
517        frames.push(MdFrame {
518            step,
519            time_fs: step as f64 * dt_fs,
520            coords: x.clone(),
521            potential_energy_kcal_mol: potential,
522            kinetic_energy_kcal_mol: ke,
523            temperature_k: temp_k,
524        });
525    }
526
527    let energy_drift_percent = if target_temp_and_tau.is_none() {
528        let e0 = frames[0].potential_energy_kcal_mol + frames[0].kinetic_energy_kcal_mol;
529        let ef = frames
530            .last()
531            .map(|f| f.potential_energy_kcal_mol + f.kinetic_energy_kcal_mol)
532            .unwrap_or(e0);
533        if e0.abs() > 1e-10 {
534            Some((ef - e0) / e0.abs() * 100.0)
535        } else {
536            None
537        }
538    } else {
539        None
540    };
541
542    Ok(MdTrajectory {
543        frames,
544        dt_fs,
545        notes: vec![format!("Velocity Verlet with {:?} backend.", backend)],
546        energy_drift_percent,
547    })
548}
549
550/// Nosé-Hoover chain thermostat for rigorous NVT sampling.
551///
552/// Implements a chain of `chain_length` thermostats coupled to velocities.
553/// Produces canonical (NVT) ensemble with correct fluctuations.
554pub fn simulate_nose_hoover(
555    smiles: &str,
556    coords: &[f64],
557    elements: &[u8],
558    n_steps: usize,
559    dt_fs: f64,
560    target_temp_k: f64,
561    thermostat_mass: f64,
562    chain_length: usize,
563    seed: u64,
564    backend: MdBackend,
565) -> Result<MdTrajectory, String> {
566    if n_steps == 0 {
567        return Err("n_steps must be > 0".to_string());
568    }
569    let n_atoms = elements.len();
570    if coords.len() != n_atoms * 3 {
571        return Err("coords length mismatch".to_string());
572    }
573
574    let masses_amu: Vec<f64> = elements.iter().map(|&z| atomic_mass_amu(z)).collect();
575    let dof = (3 * n_atoms).saturating_sub(6).max(1) as f64;
576    let target_ke = 0.5 * dof * R_GAS_KCAL_MOLK * target_temp_k;
577
578    let mut x = coords.to_vec();
579    let mut grad = vec![0.0f64; n_atoms * 3];
580    let mut potential =
581        compute_backend_energy_and_gradients(backend, smiles, elements, &x, &mut grad)?;
582
583    let mut rng = StdRng::seed_from_u64(seed);
584    let mut v = vec![0.0f64; n_atoms * 3];
585    for i in 0..n_atoms {
586        let sigma =
587            ((R_GAS_KCAL_MOLK * target_temp_k) / (masses_amu[i] * AMU_ANGFS2_TO_KCAL_MOL)).sqrt();
588        v[3 * i] = sigma * sample_standard_normal(&mut rng);
589        v[3 * i + 1] = sigma * sample_standard_normal(&mut rng);
590        v[3 * i + 2] = sigma * sample_standard_normal(&mut rng);
591    }
592
593    // Nosé-Hoover chain variables
594    let chain_len = chain_length.max(1);
595    let q = vec![thermostat_mass; chain_len]; // thermostat masses
596    let mut xi = vec![0.0f64; chain_len]; // thermostat positions
597    let mut v_xi = vec![0.0f64; chain_len]; // thermostat velocities
598
599    let (ke0, t0) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
600    let mut frames = vec![MdFrame {
601        step: 0,
602        time_fs: 0.0,
603        coords: x.clone(),
604        potential_energy_kcal_mol: potential,
605        kinetic_energy_kcal_mol: ke0,
606        temperature_k: t0,
607    }];
608
609    let dt2 = dt_fs * 0.5;
610    let dt4 = dt_fs * 0.25;
611
612    for step in 1..=n_steps {
613        let (ke, _) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
614
615        // Nosé-Hoover chain: propagate thermostat (Yoshida-Suzuki-like)
616        // Force on first thermostat
617        let mut g_xi = vec![0.0f64; chain_len];
618        g_xi[0] = (2.0 * ke - 2.0 * target_ke) / q[0];
619        for j in 1..chain_len {
620            g_xi[j] =
621                (q[j - 1] * v_xi[j - 1] * v_xi[j - 1] - R_GAS_KCAL_MOLK * target_temp_k) / q[j];
622        }
623
624        // Update thermostat velocities (outside-in)
625        if chain_len > 1 {
626            v_xi[chain_len - 1] += g_xi[chain_len - 1] * dt4;
627        }
628        for j in (0..chain_len.saturating_sub(1)).rev() {
629            let exp_factor = (-v_xi[j + 1] * dt4).exp();
630            v_xi[j] = v_xi[j] * exp_factor + g_xi[j] * dt4;
631        }
632
633        // Scale velocities
634        let scale = (-v_xi[0] * dt2).exp();
635        for vi in v.iter_mut() {
636            *vi *= scale;
637        }
638
639        // Update thermostat positions
640        for j in 0..chain_len {
641            xi[j] += v_xi[j] * dt2;
642        }
643
644        // Recompute KE after scaling
645        let (ke_scaled, _) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
646        g_xi[0] = (2.0 * ke_scaled - 2.0 * target_ke) / q[0];
647
648        // Update thermostat velocities again (inside-out)
649        for j in 0..chain_len.saturating_sub(1) {
650            let exp_factor = (-v_xi[j + 1] * dt4).exp();
651            v_xi[j] = v_xi[j] * exp_factor + g_xi[j] * dt4;
652        }
653        if chain_len > 1 {
654            // Recompute force on last thermostat
655            g_xi[chain_len - 1] = (q[chain_len - 2] * v_xi[chain_len - 2] * v_xi[chain_len - 2]
656                - R_GAS_KCAL_MOLK * target_temp_k)
657                / q[chain_len - 1];
658            v_xi[chain_len - 1] += g_xi[chain_len - 1] * dt4;
659        }
660
661        // Velocity Verlet: half-step velocity
662        for i in 0..n_atoms {
663            let m = masses_amu[i];
664            for k in 0..3 {
665                v[3 * i + k] -= 0.5 * dt_fs * grad[3 * i + k] / (m * AMU_ANGFS2_TO_KCAL_MOL);
666            }
667        }
668
669        // Position update
670        for i in 0..(n_atoms * 3) {
671            x[i] += v[i] * dt_fs;
672        }
673
674        // New forces
675        potential = compute_backend_energy_and_gradients(backend, smiles, elements, &x, &mut grad)?;
676
677        // Velocity Verlet: second half-step velocity
678        for i in 0..n_atoms {
679            let m = masses_amu[i];
680            for k in 0..3 {
681                v[3 * i + k] -= 0.5 * dt_fs * grad[3 * i + k] / (m * AMU_ANGFS2_TO_KCAL_MOL);
682            }
683        }
684
685        let (ke_final, temp_k) = kinetic_energy_and_temperature(&v, &masses_amu, n_atoms);
686
687        if !x.iter().all(|v| v.is_finite()) || !potential.is_finite() {
688            return Err(format!("NH-MD diverged at step {}", step));
689        }
690
691        frames.push(MdFrame {
692            step,
693            time_fs: step as f64 * dt_fs,
694            coords: x.clone(),
695            potential_energy_kcal_mol: potential,
696            kinetic_energy_kcal_mol: ke_final,
697            temperature_k: temp_k,
698        });
699    }
700
701    // Temperature analysis
702    let temps: Vec<f64> = frames
703        .iter()
704        .skip(frames.len() / 5)
705        .map(|f| f.temperature_k)
706        .collect();
707    let avg_temp = temps.iter().sum::<f64>() / temps.len() as f64;
708    let temp_variance =
709        temps.iter().map(|t| (t - avg_temp).powi(2)).sum::<f64>() / temps.len() as f64;
710
711    Ok(MdTrajectory {
712        frames,
713        dt_fs,
714        notes: vec![
715            format!(
716                "Nosé-Hoover chain thermostat (chain={}) with {:?} backend.",
717                chain_len, backend
718            ),
719            format!(
720                "Target T = {:.1} K, avg T = {:.1} K, σ(T) = {:.1} K",
721                target_temp_k,
722                avg_temp,
723                temp_variance.sqrt()
724            ),
725        ],
726        energy_drift_percent: None,
727    })
728}
729
730/// Export an MD trajectory to XYZ format string.
731///
732/// Each frame becomes a block with atom count, comment line (step/energy/temp),
733/// and atom coordinates.
734pub fn trajectory_to_xyz(trajectory: &MdTrajectory, elements: &[u8]) -> String {
735    let n_atoms = elements.len();
736    let mut output = String::new();
737
738    for frame in &trajectory.frames {
739        output.push_str(&format!("{}\n", n_atoms));
740        output.push_str(&format!(
741            "step={} t={:.1}fs E_pot={:.4} E_kin={:.4} T={:.1}K\n",
742            frame.step,
743            frame.time_fs,
744            frame.potential_energy_kcal_mol,
745            frame.kinetic_energy_kcal_mol,
746            frame.temperature_k,
747        ));
748        for i in 0..n_atoms {
749            let sym = element_symbol(elements[i]);
750            output.push_str(&format!(
751                "{} {:.6} {:.6} {:.6}\n",
752                sym,
753                frame.coords[3 * i],
754                frame.coords[3 * i + 1],
755                frame.coords[3 * i + 2],
756            ));
757        }
758    }
759    output
760}
761
762fn element_symbol(z: u8) -> &'static str {
763    match z {
764        1 => "H",
765        2 => "He",
766        3 => "Li",
767        4 => "Be",
768        5 => "B",
769        6 => "C",
770        7 => "N",
771        8 => "O",
772        9 => "F",
773        10 => "Ne",
774        11 => "Na",
775        12 => "Mg",
776        13 => "Al",
777        14 => "Si",
778        15 => "P",
779        16 => "S",
780        17 => "Cl",
781        18 => "Ar",
782        19 => "K",
783        20 => "Ca",
784        22 => "Ti",
785        24 => "Cr",
786        25 => "Mn",
787        26 => "Fe",
788        27 => "Co",
789        28 => "Ni",
790        29 => "Cu",
791        30 => "Zn",
792        35 => "Br",
793        44 => "Ru",
794        46 => "Pd",
795        47 => "Ag",
796        53 => "I",
797        78 => "Pt",
798        79 => "Au",
799        _ => "X",
800    }
801}