1#![allow(clippy::too_many_arguments)]
3#![allow(clippy::needless_range_loop)]
4
5pub mod alignment;
6pub mod ani;
7pub mod charges;
8pub mod clustering;
9pub mod conformer;
10pub mod dipole;
11pub mod distgeom;
12pub mod dos;
13pub mod dynamics;
14pub mod eht;
15pub mod esp;
16pub mod etkdg;
17pub mod forcefield;
18pub mod graph;
19pub mod hf;
20pub mod ir;
21pub mod materials;
22pub mod mesh;
23pub mod ml;
24pub mod nmr;
25pub mod optimization;
26pub mod pm3;
27pub mod population;
28pub mod reactivity;
29pub mod rings;
30pub mod smarts;
31pub mod smiles;
32pub mod solvation;
33pub mod stereo;
34pub mod surface;
35pub mod topology;
36pub mod transport;
37pub mod xtb;
38
39use serde::{Deserialize, Serialize};
40
41#[derive(Debug, Clone, Serialize, Deserialize)]
45pub struct ConformerResult {
46 pub smiles: String,
48 pub num_atoms: usize,
50 pub coords: Vec<f64>,
53 pub elements: Vec<u8>,
55 pub bonds: Vec<(usize, usize, String)>,
57 pub error: Option<String>,
59 pub time_ms: f64,
61}
62
63#[derive(Debug, Clone, Serialize, Deserialize)]
65pub struct ConformerConfig {
66 pub seed: u64,
68 pub num_threads: usize,
70}
71
72#[derive(Debug, Clone, Serialize, Deserialize)]
74pub struct ConformerEnsembleMember {
75 pub seed: u64,
77 pub cluster_id: Option<usize>,
79 pub coords: Vec<f64>,
81 pub energy_kcal_mol: f64,
83}
84
85#[derive(Debug, Clone, Serialize, Deserialize)]
87pub struct ConformerClusterSummary {
88 pub cluster_id: usize,
90 pub representative_seed: u64,
92 pub size: usize,
94 pub member_seeds: Vec<u64>,
96}
97
98#[derive(Debug, Clone, Serialize, Deserialize)]
100pub struct ConformerSearchResult {
101 pub generated: usize,
103 pub unique: usize,
105 pub rotatable_bonds: usize,
107 pub conformers: Vec<ConformerEnsembleMember>,
109 pub clusters: Vec<ConformerClusterSummary>,
111 pub notes: Vec<String>,
113}
114
115#[derive(Debug, Clone, Serialize, Deserialize)]
117pub struct MethodCapability {
118 pub available: bool,
120 pub confidence: eht::SupportLevel,
122 pub unsupported_elements: Vec<u8>,
124 pub warnings: Vec<String>,
126}
127
128#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
130#[serde(rename_all = "snake_case")]
131pub enum ScientificMethod {
132 Embed,
133 Uff,
134 Eht,
135 Pm3,
136 Xtb,
137 Mmff94,
138 Ani,
139 Hf3c,
140}
141
142#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
144#[serde(rename_all = "snake_case")]
145pub enum PropertyRequest {
146 Geometry,
147 ForceFieldEnergy,
148 Orbitals,
149 Population,
150 OrbitalGrid,
151}
152
153#[derive(Debug, Clone, Serialize, Deserialize)]
155pub struct MethodMetadata {
156 pub method: ScientificMethod,
157 pub available: bool,
158 pub confidence: eht::SupportLevel,
159 pub confidence_score: f64,
160 pub limitations: Vec<String>,
161 pub warnings: Vec<String>,
162}
163
164#[derive(Debug, Clone, Serialize, Deserialize)]
166pub struct PropertyMethodPlan {
167 pub property: PropertyRequest,
168 pub recommended: Option<ScientificMethod>,
169 pub fallback: Option<ScientificMethod>,
170 pub rationale: Vec<String>,
171 pub methods: Vec<MethodMetadata>,
172}
173
174#[derive(Debug, Clone, Serialize, Deserialize)]
176pub struct SystemMethodPlan {
177 pub capabilities: SystemCapabilities,
178 pub geometry: PropertyMethodPlan,
179 pub force_field_energy: PropertyMethodPlan,
180 pub orbitals: PropertyMethodPlan,
181 pub population: PropertyMethodPlan,
182 pub orbital_grid: PropertyMethodPlan,
183}
184
185#[derive(Debug, Clone, Serialize, Deserialize)]
187pub struct SystemCapabilities {
188 pub embed: MethodCapability,
189 pub uff: MethodCapability,
190 pub eht: MethodCapability,
191 pub population: MethodCapability,
192 pub orbital_grid: MethodCapability,
193}
194
195#[derive(Debug, Clone, Serialize, Deserialize)]
197#[serde(tag = "mode", rename_all = "snake_case")]
198pub enum ElectronicWorkflowResult {
199 Eht {
200 result: eht::EhtResult,
201 },
202 UffFallback {
203 energy_kcal_mol: f64,
204 reason: String,
205 support: eht::EhtSupport,
206 },
207}
208
209#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
211#[serde(rename_all = "snake_case")]
212pub enum MethodComparisonStatus {
213 Success,
214 Unavailable,
215 Error,
216}
217
218#[derive(Debug, Clone, Serialize, Deserialize)]
220#[serde(tag = "kind", rename_all = "snake_case")]
221pub enum MethodComparisonPayload {
222 Eht {
223 homo_energy: f64,
224 lumo_energy: f64,
225 gap: f64,
226 support: eht::EhtSupport,
227 },
228 Uff {
229 energy_kcal_mol: f64,
230 },
231}
232
233#[derive(Debug, Clone, Serialize, Deserialize)]
235pub struct MethodComparisonEntry {
236 pub method: ScientificMethod,
237 pub status: MethodComparisonStatus,
238 pub available: bool,
239 pub confidence: eht::SupportLevel,
240 pub confidence_score: f64,
241 pub warnings: Vec<String>,
242 pub limitations: Vec<String>,
243 pub payload: Option<MethodComparisonPayload>,
244 pub error: Option<String>,
245}
246
247#[derive(Debug, Clone, Serialize, Deserialize)]
249pub struct MethodComparisonResult {
250 pub plan: SystemMethodPlan,
251 pub comparisons: Vec<MethodComparisonEntry>,
252}
253
254#[derive(Debug, Clone, Serialize, Deserialize)]
256pub struct AromaticityAnalysis {
257 pub aromatic_atoms: Vec<bool>,
258 pub aromatic_bonds: Vec<(usize, usize)>,
259}
260
261#[derive(Debug, Clone, Serialize, Deserialize)]
263pub struct StereocenterAnalysis {
264 pub tagged_tetrahedral_centers: Vec<usize>,
265 pub inferred_tetrahedral_centers: Vec<usize>,
266}
267
268#[derive(Debug, Clone, Serialize, Deserialize)]
270pub struct GraphFeatureAnalysis {
271 pub aromaticity: AromaticityAnalysis,
272 pub stereocenters: StereocenterAnalysis,
273}
274
275impl Default for ConformerConfig {
276 fn default() -> Self {
277 Self {
278 seed: 42,
279 num_threads: 0,
280 }
281 }
282}
283
284pub fn version() -> String {
288 format!("sci-form {}", env!("CARGO_PKG_VERSION"))
289}
290
291pub fn get_eht_support(elements: &[u8]) -> eht::EhtSupport {
293 eht::analyze_eht_support(elements)
294}
295
296fn is_uff_element_supported(z: u8) -> bool {
297 matches!(
298 z,
299 1 | 5
300 | 6
301 | 7
302 | 8
303 | 9
304 | 14
305 | 15
306 | 16
307 | 17
308 | 22
309 | 23
310 | 24
311 | 25
312 | 26
313 | 27
314 | 28
315 | 29
316 | 30
317 | 32
318 | 33
319 | 34
320 | 35
321 | 42
322 | 46
323 | 47
324 | 50
325 | 51
326 | 52
327 | 53
328 | 78
329 | 79
330 )
331}
332
333fn unique_sorted_unsupported(elements: &[u8], pred: impl Fn(u8) -> bool) -> Vec<u8> {
334 let mut out: Vec<u8> = elements.iter().copied().filter(|&z| !pred(z)).collect();
335 out.sort_unstable();
336 out.dedup();
337 out
338}
339
340pub fn get_system_capabilities(elements: &[u8]) -> SystemCapabilities {
342 let eht_support = get_eht_support(elements);
343 let uff_unsupported = unique_sorted_unsupported(elements, is_uff_element_supported);
344
345 let embed = MethodCapability {
346 available: !elements.is_empty(),
347 confidence: eht::SupportLevel::Experimental,
348 unsupported_elements: Vec::new(),
349 warnings: vec![
350 "Embed capability is inferred from element presence only; final success still depends on full molecular graph and geometry constraints.".to_string(),
351 ],
352 };
353
354 let uff = if uff_unsupported.is_empty() {
355 MethodCapability {
356 available: true,
357 confidence: eht::SupportLevel::Supported,
358 unsupported_elements: Vec::new(),
359 warnings: Vec::new(),
360 }
361 } else {
362 MethodCapability {
363 available: false,
364 confidence: eht::SupportLevel::Unsupported,
365 unsupported_elements: uff_unsupported.clone(),
366 warnings: vec![format!(
367 "UFF atom typing is unavailable for elements {:?}.",
368 uff_unsupported
369 )],
370 }
371 };
372
373 let eht = MethodCapability {
374 available: eht_support.unsupported_elements.is_empty(),
375 confidence: eht_support.level,
376 unsupported_elements: eht_support.unsupported_elements.clone(),
377 warnings: eht_support.warnings.clone(),
378 };
379
380 let population = MethodCapability {
381 available: eht.available,
382 confidence: eht.confidence,
383 unsupported_elements: eht.unsupported_elements.clone(),
384 warnings: eht.warnings.clone(),
385 };
386
387 let orbital_grid = MethodCapability {
388 available: eht.available,
389 confidence: eht.confidence,
390 unsupported_elements: eht.unsupported_elements.clone(),
391 warnings: eht.warnings.clone(),
392 };
393
394 SystemCapabilities {
395 embed,
396 uff,
397 eht,
398 population,
399 orbital_grid,
400 }
401}
402
403fn confidence_score_for_method(method: ScientificMethod, capability: &MethodCapability) -> f64 {
404 if !capability.available {
405 return 0.0;
406 }
407
408 match method {
409 ScientificMethod::Embed => 0.8,
410 ScientificMethod::Uff | ScientificMethod::Mmff94 => match capability.confidence {
411 eht::SupportLevel::Supported => 0.95,
412 eht::SupportLevel::Experimental => 0.75,
413 eht::SupportLevel::Unsupported => 0.0,
414 },
415 ScientificMethod::Eht | ScientificMethod::Pm3 | ScientificMethod::Xtb => {
416 match capability.confidence {
417 eht::SupportLevel::Supported => 0.95,
418 eht::SupportLevel::Experimental => 0.6,
419 eht::SupportLevel::Unsupported => 0.0,
420 }
421 }
422 ScientificMethod::Ani => match capability.confidence {
423 eht::SupportLevel::Supported => 0.90,
424 eht::SupportLevel::Experimental => 0.7,
425 eht::SupportLevel::Unsupported => 0.0,
426 },
427 ScientificMethod::Hf3c => match capability.confidence {
428 eht::SupportLevel::Supported => 0.85,
429 eht::SupportLevel::Experimental => 0.65,
430 eht::SupportLevel::Unsupported => 0.0,
431 },
432 }
433}
434
435fn build_method_metadata(
436 method: ScientificMethod,
437 capability: &MethodCapability,
438 extra_limitations: &[&str],
439) -> MethodMetadata {
440 let mut limitations: Vec<String> = extra_limitations.iter().map(|s| s.to_string()).collect();
441
442 if !capability.unsupported_elements.is_empty() {
443 limitations.push(format!(
444 "Unsupported elements for this method: {:?}.",
445 capability.unsupported_elements
446 ));
447 }
448
449 if matches!(method, ScientificMethod::Eht)
450 && matches!(capability.confidence, eht::SupportLevel::Experimental)
451 {
452 limitations.push(
453 "Transition-metal EHT parameters remain provisional and should be treated as experimental."
454 .to_string(),
455 );
456 }
457
458 MethodMetadata {
459 method,
460 available: capability.available,
461 confidence: capability.confidence,
462 confidence_score: confidence_score_for_method(method, capability),
463 limitations,
464 warnings: capability.warnings.clone(),
465 }
466}
467
468fn build_property_plan(
469 property: PropertyRequest,
470 recommended: Option<ScientificMethod>,
471 fallback: Option<ScientificMethod>,
472 rationale: Vec<String>,
473 methods: Vec<MethodMetadata>,
474) -> PropertyMethodPlan {
475 PropertyMethodPlan {
476 property,
477 recommended,
478 fallback,
479 rationale,
480 methods,
481 }
482}
483
484pub fn get_system_method_plan(elements: &[u8]) -> SystemMethodPlan {
486 let capabilities = get_system_capabilities(elements);
487
488 let geometry_method = build_method_metadata(
489 ScientificMethod::Embed,
490 &capabilities.embed,
491 &["Geometry generation still depends on graph topology, stereochemistry, and embedding constraints."],
492 );
493 let geometry = build_property_plan(
494 PropertyRequest::Geometry,
495 geometry_method.available.then_some(ScientificMethod::Embed),
496 None,
497 vec!["Embedding is the top-level geometry generation path in sci-form.".to_string()],
498 vec![geometry_method],
499 );
500
501 let uff_method = build_method_metadata(
502 ScientificMethod::Uff,
503 &capabilities.uff,
504 &["This recommendation applies to force-field energy evaluation, not molecular orbital analysis."],
505 );
506 let force_field_energy = build_property_plan(
507 PropertyRequest::ForceFieldEnergy,
508 uff_method.available.then_some(ScientificMethod::Uff),
509 None,
510 vec![
511 "UFF is the top-level force-field energy path when atom typing is available."
512 .to_string(),
513 ],
514 vec![uff_method],
515 );
516
517 let eht_method = build_method_metadata(
518 ScientificMethod::Eht,
519 &capabilities.eht,
520 &["EHT is the only current top-level orbital method in sci-form."],
521 );
522 let orbitals = build_property_plan(
523 PropertyRequest::Orbitals,
524 eht_method.available.then_some(ScientificMethod::Eht),
525 None,
526 vec!["Orbital energies and MO coefficients are produced by the EHT workflow.".to_string()],
527 vec![eht_method.clone()],
528 );
529
530 let population_method = build_method_metadata(
531 ScientificMethod::Eht,
532 &capabilities.population,
533 &["Population analysis is derived from the EHT density and overlap matrices."],
534 );
535 let population = build_property_plan(
536 PropertyRequest::Population,
537 population_method.available.then_some(ScientificMethod::Eht),
538 None,
539 vec!["Population analysis currently requires a successful EHT calculation.".to_string()],
540 vec![population_method],
541 );
542
543 let orbital_grid_method = build_method_metadata(
544 ScientificMethod::Eht,
545 &capabilities.orbital_grid,
546 &["Orbital-grid rendering currently depends on EHT molecular-orbital coefficients."],
547 );
548 let orbital_grid = build_property_plan(
549 PropertyRequest::OrbitalGrid,
550 orbital_grid_method
551 .available
552 .then_some(ScientificMethod::Eht),
553 None,
554 vec![
555 "Orbital-grid generation currently requires a successful EHT calculation.".to_string(),
556 ],
557 vec![orbital_grid_method],
558 );
559
560 SystemMethodPlan {
561 capabilities,
562 geometry,
563 force_field_energy,
564 orbitals,
565 population,
566 orbital_grid,
567 }
568}
569
570pub fn embed(smiles: &str, seed: u64) -> ConformerResult {
572 #[cfg(not(target_arch = "wasm32"))]
573 let start = std::time::Instant::now();
574
575 let mol = match graph::Molecule::from_smiles(smiles) {
576 Ok(m) => m,
577 Err(e) => {
578 return ConformerResult {
579 smiles: smiles.to_string(),
580 num_atoms: 0,
581 coords: vec![],
582 elements: vec![],
583 bonds: vec![],
584 error: Some(e),
585 #[cfg(not(target_arch = "wasm32"))]
586 time_ms: start.elapsed().as_secs_f64() * 1000.0,
587 #[cfg(target_arch = "wasm32")]
588 time_ms: 0.0,
589 };
590 }
591 };
592
593 let n = mol.graph.node_count();
594 let elements: Vec<u8> = (0..n)
595 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
596 .collect();
597 let bonds: Vec<(usize, usize, String)> = mol
598 .graph
599 .edge_indices()
600 .map(|e| {
601 let (a, b) = mol.graph.edge_endpoints(e).unwrap();
602 let order = match mol.graph[e].order {
603 graph::BondOrder::Single => "SINGLE",
604 graph::BondOrder::Double => "DOUBLE",
605 graph::BondOrder::Triple => "TRIPLE",
606 graph::BondOrder::Aromatic => "AROMATIC",
607 graph::BondOrder::Unknown => "UNKNOWN",
608 };
609 (a.index(), b.index(), order.to_string())
610 })
611 .collect();
612
613 let max_seed_retries = 3u64;
616 let mut last_err = String::new();
617 for retry in 0..max_seed_retries {
618 let current_seed = seed.wrapping_add(retry.wrapping_mul(997));
619 match conformer::generate_3d_conformer(&mol, current_seed) {
620 Ok(coords) => {
621 let mut flat = Vec::with_capacity(n * 3);
622 for i in 0..n {
623 flat.push(coords[(i, 0)] as f64);
624 flat.push(coords[(i, 1)] as f64);
625 flat.push(coords[(i, 2)] as f64);
626 }
627 return ConformerResult {
628 smiles: smiles.to_string(),
629 num_atoms: n,
630 coords: flat,
631 elements,
632 bonds,
633 error: None,
634 #[cfg(not(target_arch = "wasm32"))]
635 time_ms: start.elapsed().as_secs_f64() * 1000.0,
636 #[cfg(target_arch = "wasm32")]
637 time_ms: 0.0,
638 };
639 }
640 Err(e) => {
641 last_err = e;
642 }
643 }
644 }
645 ConformerResult {
646 smiles: smiles.to_string(),
647 num_atoms: n,
648 coords: vec![],
649 elements,
650 bonds,
651 error: Some(last_err),
652 #[cfg(not(target_arch = "wasm32"))]
653 time_ms: start.elapsed().as_secs_f64() * 1000.0,
654 #[cfg(target_arch = "wasm32")]
655 time_ms: 0.0,
656 }
657}
658
659#[cfg(feature = "parallel")]
664pub fn embed_batch(smiles_list: &[&str], config: &ConformerConfig) -> Vec<ConformerResult> {
665 use rayon::prelude::*;
666
667 if config.num_threads > 0 {
668 let pool = rayon::ThreadPoolBuilder::new()
669 .num_threads(config.num_threads)
670 .build()
671 .unwrap();
672 pool.install(|| {
673 smiles_list
674 .par_iter()
675 .map(|smi| embed(smi, config.seed))
676 .collect()
677 })
678 } else {
679 smiles_list
680 .par_iter()
681 .map(|smi| embed(smi, config.seed))
682 .collect()
683 }
684}
685
686#[cfg(not(feature = "parallel"))]
688pub fn embed_batch(smiles_list: &[&str], config: &ConformerConfig) -> Vec<ConformerResult> {
689 smiles_list
690 .iter()
691 .map(|smi| embed(smi, config.seed))
692 .collect()
693}
694
695pub fn parse(smiles: &str) -> Result<graph::Molecule, String> {
697 graph::Molecule::from_smiles(smiles)
698}
699
700pub fn compute_charges(smiles: &str) -> Result<charges::gasteiger::ChargeResult, String> {
705 let mol = graph::Molecule::from_smiles(smiles)?;
706 let n = mol.graph.node_count();
707 let elements: Vec<u8> = (0..n)
708 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
709 .collect();
710 let formal_charges: Vec<i8> = (0..n)
711 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].formal_charge)
712 .collect();
713 let bonds: Vec<(usize, usize)> = mol
714 .graph
715 .edge_indices()
716 .map(|e| {
717 let (a, b) = mol.graph.edge_endpoints(e).unwrap();
718 (a.index(), b.index())
719 })
720 .collect();
721 charges::gasteiger::gasteiger_marsili_charges(&elements, &bonds, &formal_charges, 6)
722}
723
724pub fn compute_sasa(
728 elements: &[u8],
729 coords_flat: &[f64],
730 probe_radius: Option<f64>,
731) -> Result<surface::sasa::SasaResult, String> {
732 if coords_flat.len() != elements.len() * 3 {
733 return Err(format!(
734 "coords length {} != 3 * elements {}",
735 coords_flat.len(),
736 elements.len()
737 ));
738 }
739 let positions: Vec<[f64; 3]> = coords_flat
740 .chunks_exact(3)
741 .map(|c| [c[0], c[1], c[2]])
742 .collect();
743
744 #[cfg(feature = "parallel")]
745 {
746 Ok(surface::sasa::compute_sasa_parallel(
747 elements,
748 &positions,
749 probe_radius,
750 None,
751 ))
752 }
753
754 #[cfg(not(feature = "parallel"))]
755 {
756 Ok(surface::sasa::compute_sasa(
757 elements,
758 &positions,
759 probe_radius,
760 None,
761 ))
762 }
763}
764
765pub fn compute_population(
767 elements: &[u8],
768 positions: &[[f64; 3]],
769) -> Result<population::PopulationResult, String> {
770 let eht_result = eht::solve_eht(elements, positions, None)?;
771 Ok(population::compute_population(
772 elements,
773 positions,
774 &eht_result.coefficients,
775 eht_result.n_electrons,
776 ))
777}
778
779pub fn compute_dipole(
781 elements: &[u8],
782 positions: &[[f64; 3]],
783) -> Result<dipole::DipoleResult, String> {
784 let eht_result = eht::solve_eht(elements, positions, None)?;
785 Ok(dipole::compute_dipole_from_eht(
786 elements,
787 positions,
788 &eht_result.coefficients,
789 eht_result.n_electrons,
790 ))
791}
792
793pub fn compute_frontier_descriptors(
795 elements: &[u8],
796 positions: &[[f64; 3]],
797) -> Result<reactivity::FrontierDescriptors, String> {
798 let eht_result = eht::solve_eht(elements, positions, None)?;
799 Ok(reactivity::compute_frontier_descriptors(
800 elements,
801 positions,
802 &eht_result,
803 ))
804}
805
806pub fn compute_fukui_descriptors(
808 elements: &[u8],
809 positions: &[[f64; 3]],
810) -> Result<reactivity::FukuiDescriptors, String> {
811 let eht_result = eht::solve_eht(elements, positions, None)?;
812 Ok(reactivity::compute_fukui_descriptors(
813 elements,
814 positions,
815 &eht_result,
816 ))
817}
818
819pub fn compute_reactivity_ranking(
821 elements: &[u8],
822 positions: &[[f64; 3]],
823) -> Result<reactivity::ReactivityRanking, String> {
824 let eht_result = eht::solve_eht(elements, positions, None)?;
825 let fukui = reactivity::compute_fukui_descriptors(elements, positions, &eht_result);
826 let pop = population::compute_population(
827 elements,
828 positions,
829 &eht_result.coefficients,
830 eht_result.n_electrons,
831 );
832 Ok(reactivity::rank_reactivity_sites(
833 &fukui,
834 &pop.mulliken_charges,
835 ))
836}
837
838pub fn compute_uv_vis_spectrum(
840 elements: &[u8],
841 positions: &[[f64; 3]],
842 sigma: f64,
843 e_min: f64,
844 e_max: f64,
845 n_points: usize,
846) -> Result<reactivity::UvVisSpectrum, String> {
847 let eht_result = eht::solve_eht(elements, positions, None)?;
848 Ok(reactivity::compute_uv_vis_like_spectrum(
849 &eht_result,
850 sigma,
851 e_min,
852 e_max,
853 n_points,
854 ))
855}
856
857pub fn compute_bond_orders(
859 elements: &[u8],
860 positions: &[[f64; 3]],
861) -> Result<population::BondOrderResult, String> {
862 let eht_result = eht::solve_eht(elements, positions, None)?;
863 Ok(population::compute_bond_orders(
864 elements,
865 positions,
866 &eht_result.coefficients,
867 eht_result.n_electrons,
868 ))
869}
870
871pub fn compute_topology(
873 elements: &[u8],
874 positions: &[[f64; 3]],
875) -> topology::TopologyAnalysisResult {
876 topology::analyze_topology(elements, positions)
877}
878
879pub fn analyze_graph_features(smiles: &str) -> Result<GraphFeatureAnalysis, String> {
881 use petgraph::visit::EdgeRef;
882
883 let mol = parse(smiles)?;
884 let n_atoms = mol.graph.node_count();
885 let mut aromatic_atoms = vec![false; n_atoms];
886 let mut aromatic_bonds = Vec::new();
887
888 for edge in mol.graph.edge_references() {
889 if matches!(edge.weight().order, graph::BondOrder::Aromatic) {
890 let i = edge.source().index();
891 let j = edge.target().index();
892 aromatic_atoms[i] = true;
893 aromatic_atoms[j] = true;
894 aromatic_bonds.push((i, j));
895 }
896 }
897
898 let mut tagged_tetrahedral_centers = Vec::new();
899 let mut inferred_tetrahedral_centers = Vec::new();
900 for i in 0..n_atoms {
901 let idx = petgraph::graph::NodeIndex::new(i);
902 let atom = &mol.graph[idx];
903 if matches!(
904 atom.chiral_tag,
905 graph::ChiralType::TetrahedralCW | graph::ChiralType::TetrahedralCCW
906 ) {
907 tagged_tetrahedral_centers.push(i);
908 }
909
910 let neighbors: Vec<_> = mol.graph.neighbors(idx).collect();
911 if neighbors.len() == 4 && matches!(atom.hybridization, graph::Hybridization::SP3) {
912 let mut signature: Vec<u8> = neighbors.iter().map(|n| mol.graph[*n].element).collect();
913 signature.sort_unstable();
914 signature.dedup();
915 if signature.len() >= 3 {
916 inferred_tetrahedral_centers.push(i);
917 }
918 }
919 }
920
921 Ok(GraphFeatureAnalysis {
922 aromaticity: AromaticityAnalysis {
923 aromatic_atoms,
924 aromatic_bonds,
925 },
926 stereocenters: StereocenterAnalysis {
927 tagged_tetrahedral_centers,
928 inferred_tetrahedral_centers,
929 },
930 })
931}
932
933pub fn compute_eht_or_uff_fallback(
938 smiles: &str,
939 elements: &[u8],
940 positions: &[[f64; 3]],
941 allow_experimental_eht: bool,
942) -> Result<ElectronicWorkflowResult, String> {
943 let support = get_eht_support(elements);
944 let should_fallback = match support.level {
945 eht::SupportLevel::Unsupported => true,
946 eht::SupportLevel::Experimental => !allow_experimental_eht,
947 eht::SupportLevel::Supported => false,
948 };
949
950 if should_fallback {
951 let coords_flat: Vec<f64> = positions.iter().flat_map(|p| p.iter().copied()).collect();
952 let energy = compute_uff_energy(smiles, &coords_flat).map_err(|e| {
953 format!(
954 "EHT is not appropriate for this system (support: {:?}) and UFF fallback failed: {}",
955 support.level, e
956 )
957 })?;
958 return Ok(ElectronicWorkflowResult::UffFallback {
959 energy_kcal_mol: energy,
960 reason: if matches!(support.level, eht::SupportLevel::Unsupported) {
961 "EHT unsupported for one or more elements; routed to UFF-only workflow.".to_string()
962 } else {
963 "EHT confidence is experimental and experimental mode is disabled; routed to UFF-only workflow."
964 .to_string()
965 },
966 support,
967 });
968 }
969
970 let result = eht::solve_eht(elements, positions, None)?;
971 Ok(ElectronicWorkflowResult::Eht { result })
972}
973
974pub fn compare_methods(
979 smiles: &str,
980 elements: &[u8],
981 positions: &[[f64; 3]],
982 allow_experimental_eht: bool,
983) -> MethodComparisonResult {
984 let plan = get_system_method_plan(elements);
985 let mut comparisons = Vec::new();
986
987 let coords_flat: Vec<f64> = positions.iter().flat_map(|p| p.iter().copied()).collect();
988
989 {
990 let meta = build_method_metadata(
991 ScientificMethod::Uff,
992 &plan.capabilities.uff,
993 &["Comparison uses UFF force-field energy as the UFF observable."],
994 );
995 if !meta.available {
996 comparisons.push(MethodComparisonEntry {
997 method: ScientificMethod::Uff,
998 status: MethodComparisonStatus::Unavailable,
999 available: false,
1000 confidence: meta.confidence,
1001 confidence_score: meta.confidence_score,
1002 warnings: meta.warnings,
1003 limitations: meta.limitations,
1004 payload: None,
1005 error: Some("UFF is unavailable for this element set.".to_string()),
1006 });
1007 } else {
1008 match compute_uff_energy(smiles, &coords_flat) {
1009 Ok(energy) => comparisons.push(MethodComparisonEntry {
1010 method: ScientificMethod::Uff,
1011 status: MethodComparisonStatus::Success,
1012 available: true,
1013 confidence: meta.confidence,
1014 confidence_score: meta.confidence_score,
1015 warnings: meta.warnings,
1016 limitations: meta.limitations,
1017 payload: Some(MethodComparisonPayload::Uff {
1018 energy_kcal_mol: energy,
1019 }),
1020 error: None,
1021 }),
1022 Err(err) => comparisons.push(MethodComparisonEntry {
1023 method: ScientificMethod::Uff,
1024 status: MethodComparisonStatus::Error,
1025 available: true,
1026 confidence: meta.confidence,
1027 confidence_score: meta.confidence_score,
1028 warnings: meta.warnings,
1029 limitations: meta.limitations,
1030 payload: None,
1031 error: Some(err),
1032 }),
1033 }
1034 }
1035 }
1036
1037 {
1038 let meta = build_method_metadata(
1039 ScientificMethod::Eht,
1040 &plan.capabilities.eht,
1041 &["Comparison uses frontier orbital energies and gap as the EHT observable."],
1042 );
1043
1044 if !meta.available {
1045 comparisons.push(MethodComparisonEntry {
1046 method: ScientificMethod::Eht,
1047 status: MethodComparisonStatus::Unavailable,
1048 available: false,
1049 confidence: meta.confidence,
1050 confidence_score: meta.confidence_score,
1051 warnings: meta.warnings,
1052 limitations: meta.limitations,
1053 payload: None,
1054 error: Some("EHT is unavailable for this element set.".to_string()),
1055 });
1056 } else if matches!(meta.confidence, eht::SupportLevel::Experimental)
1057 && !allow_experimental_eht
1058 {
1059 comparisons.push(MethodComparisonEntry {
1060 method: ScientificMethod::Eht,
1061 status: MethodComparisonStatus::Unavailable,
1062 available: true,
1063 confidence: meta.confidence,
1064 confidence_score: meta.confidence_score,
1065 warnings: meta.warnings,
1066 limitations: meta.limitations,
1067 payload: None,
1068 error: Some(
1069 "EHT confidence is experimental and allow_experimental_eht=false.".to_string(),
1070 ),
1071 });
1072 } else {
1073 match eht::solve_eht(elements, positions, None) {
1074 Ok(result) => comparisons.push(MethodComparisonEntry {
1075 method: ScientificMethod::Eht,
1076 status: MethodComparisonStatus::Success,
1077 available: true,
1078 confidence: meta.confidence,
1079 confidence_score: meta.confidence_score,
1080 warnings: meta.warnings,
1081 limitations: meta.limitations,
1082 payload: Some(MethodComparisonPayload::Eht {
1083 homo_energy: result.homo_energy,
1084 lumo_energy: result.lumo_energy,
1085 gap: result.gap,
1086 support: result.support,
1087 }),
1088 error: None,
1089 }),
1090 Err(err) => comparisons.push(MethodComparisonEntry {
1091 method: ScientificMethod::Eht,
1092 status: MethodComparisonStatus::Error,
1093 available: true,
1094 confidence: meta.confidence,
1095 confidence_score: meta.confidence_score,
1096 warnings: meta.warnings,
1097 limitations: meta.limitations,
1098 payload: None,
1099 error: Some(err),
1100 }),
1101 }
1102 }
1103 }
1104
1105 MethodComparisonResult { plan, comparisons }
1106}
1107
1108pub fn compute_esp(
1110 elements: &[u8],
1111 positions: &[[f64; 3]],
1112 spacing: f64,
1113 padding: f64,
1114) -> Result<esp::EspGrid, String> {
1115 let pop = compute_population(elements, positions)?;
1116 #[cfg(feature = "parallel")]
1117 {
1118 Ok(esp::compute_esp_grid_parallel(
1119 elements,
1120 positions,
1121 &pop.mulliken_charges,
1122 spacing,
1123 padding,
1124 ))
1125 }
1126
1127 #[cfg(not(feature = "parallel"))]
1128 {
1129 Ok(esp::compute_esp_grid(
1130 elements,
1131 positions,
1132 &pop.mulliken_charges,
1133 spacing,
1134 padding,
1135 ))
1136 }
1137}
1138
1139pub fn compute_dos(
1141 elements: &[u8],
1142 positions: &[[f64; 3]],
1143 sigma: f64,
1144 e_min: f64,
1145 e_max: f64,
1146 n_points: usize,
1147) -> Result<dos::DosResult, String> {
1148 let eht_result = eht::solve_eht(elements, positions, None)?;
1149 let flat: Vec<f64> = positions.iter().flat_map(|p| p.iter().copied()).collect();
1150
1151 #[cfg(feature = "parallel")]
1152 {
1153 Ok(dos::compute_pdos_parallel(
1154 elements,
1155 &flat,
1156 &eht_result.energies,
1157 &eht_result.coefficients,
1158 eht_result.n_electrons,
1159 sigma,
1160 e_min,
1161 e_max,
1162 n_points,
1163 ))
1164 }
1165
1166 #[cfg(not(feature = "parallel"))]
1167 {
1168 Ok(dos::compute_pdos(
1169 elements,
1170 &flat,
1171 &eht_result.energies,
1172 &eht_result.coefficients,
1173 eht_result.n_electrons,
1174 sigma,
1175 e_min,
1176 e_max,
1177 n_points,
1178 ))
1179 }
1180}
1181
1182pub fn compute_rmsd(coords: &[f64], reference: &[f64]) -> f64 {
1184 alignment::compute_rmsd(coords, reference)
1185}
1186
1187pub fn compute_uff_energy(smiles: &str, coords: &[f64]) -> Result<f64, String> {
1189 let mol = graph::Molecule::from_smiles(smiles)?;
1190 let n = mol.graph.node_count();
1191 if coords.len() != n * 3 {
1192 return Err(format!("coords length {} != 3 * atoms {}", coords.len(), n));
1193 }
1194 let ff = forcefield::builder::build_uff_force_field(&mol);
1195 let mut gradient = vec![0.0f64; n * 3];
1196 let energy = ff.compute_system_energy_and_gradients(coords, &mut gradient);
1197 Ok(energy)
1198}
1199
1200pub fn compute_mmff94_energy(smiles: &str, coords: &[f64]) -> Result<f64, String> {
1207 let mol = graph::Molecule::from_smiles(smiles)?;
1208 let n = mol.graph.node_count();
1209 if coords.len() != n * 3 {
1210 return Err(format!("coords length {} != 3 * atoms {}", coords.len(), n));
1211 }
1212 let elements: Vec<u8> = (0..n)
1213 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
1214 .collect();
1215 let bonds: Vec<(usize, usize, u8)> = mol
1216 .graph
1217 .edge_indices()
1218 .map(|e| {
1219 let (a, b) = mol.graph.edge_endpoints(e).unwrap();
1220 let order = match mol.graph[e].order {
1221 graph::BondOrder::Single => 1u8,
1222 graph::BondOrder::Double => 2,
1223 graph::BondOrder::Triple => 3,
1224 graph::BondOrder::Aromatic => 2,
1225 graph::BondOrder::Unknown => 1,
1226 };
1227 (a.index(), b.index(), order)
1228 })
1229 .collect();
1230 let terms = forcefield::mmff94::Mmff94Builder::build(&elements, &bonds);
1231 let (energy, _grad) = forcefield::mmff94::Mmff94Builder::total_energy(&terms, coords);
1232 Ok(energy)
1233}
1234
1235pub fn compute_pm3(elements: &[u8], positions: &[[f64; 3]]) -> Result<pm3::Pm3Result, String> {
1242 pm3::solve_pm3(elements, positions)
1243}
1244
1245pub fn compute_xtb(elements: &[u8], positions: &[[f64; 3]]) -> Result<xtb::XtbResult, String> {
1252 xtb::solve_xtb(elements, positions)
1253}
1254
1255pub fn compute_stda_uvvis(
1262 elements: &[u8],
1263 positions: &[[f64; 3]],
1264 sigma: f64,
1265 e_min: f64,
1266 e_max: f64,
1267 n_points: usize,
1268 broadening: reactivity::BroadeningType,
1269) -> Result<reactivity::StdaUvVisSpectrum, String> {
1270 reactivity::compute_stda_uvvis_spectrum(
1271 elements, positions, sigma, e_min, e_max, n_points, broadening,
1272 )
1273}
1274
1275pub fn compute_vibrational_analysis(
1282 elements: &[u8],
1283 positions: &[[f64; 3]],
1284 method: &str,
1285 step_size: Option<f64>,
1286) -> Result<ir::VibrationalAnalysis, String> {
1287 let hessian_method = match method {
1288 "eht" => ir::HessianMethod::Eht,
1289 "pm3" => ir::HessianMethod::Pm3,
1290 "xtb" => ir::HessianMethod::Xtb,
1291 _ => return Err(format!("Unknown method '{}', use eht/pm3/xtb", method)),
1292 };
1293 ir::compute_vibrational_analysis(elements, positions, hessian_method, step_size)
1294}
1295
1296pub fn compute_ir_spectrum(
1302 analysis: &ir::VibrationalAnalysis,
1303 gamma: f64,
1304 wn_min: f64,
1305 wn_max: f64,
1306 n_points: usize,
1307) -> ir::IrSpectrum {
1308 ir::compute_ir_spectrum(analysis, gamma, wn_min, wn_max, n_points)
1309}
1310
1311pub fn predict_nmr_shifts(smiles: &str) -> Result<nmr::NmrShiftResult, String> {
1313 let mol = graph::Molecule::from_smiles(smiles)?;
1314 Ok(nmr::predict_chemical_shifts(&mol))
1315}
1316
1317pub fn predict_nmr_couplings(
1322 smiles: &str,
1323 positions: &[[f64; 3]],
1324) -> Result<Vec<nmr::JCoupling>, String> {
1325 let mol = graph::Molecule::from_smiles(smiles)?;
1326 Ok(nmr::predict_j_couplings(&mol, positions))
1327}
1328
1329pub fn compute_nmr_spectrum(
1336 smiles: &str,
1337 nucleus: &str,
1338 gamma: f64,
1339 ppm_min: f64,
1340 ppm_max: f64,
1341 n_points: usize,
1342) -> Result<nmr::NmrSpectrum, String> {
1343 let mol = graph::Molecule::from_smiles(smiles)?;
1344 let shifts = nmr::predict_chemical_shifts(&mol);
1345 let couplings = nmr::predict_j_couplings(&mol, &[]);
1346 let nuc = match nucleus {
1347 "1H" | "H1" | "h1" | "1h" | "proton" => nmr::NmrNucleus::H1,
1348 "13C" | "C13" | "c13" | "13c" | "carbon" => nmr::NmrNucleus::C13,
1349 "19F" | "F19" | "f19" | "19f" | "fluorine" => nmr::NmrNucleus::F19,
1350 "31P" | "P31" | "p31" | "31p" | "phosphorus" => nmr::NmrNucleus::P31,
1351 "15N" | "N15" | "n15" | "15n" | "nitrogen" => nmr::NmrNucleus::N15,
1352 "11B" | "B11" | "b11" | "11b" | "boron" => nmr::NmrNucleus::B11,
1353 "29Si" | "Si29" | "si29" | "29si" | "silicon" => nmr::NmrNucleus::Si29,
1354 "77Se" | "Se77" | "se77" | "77se" | "selenium" => nmr::NmrNucleus::Se77,
1355 "17O" | "O17" | "o17" | "17o" | "oxygen" => nmr::NmrNucleus::O17,
1356 "33S" | "S33" | "s33" | "33s" | "sulfur" => nmr::NmrNucleus::S33,
1357 _ => {
1358 return Err(format!(
1359 "Unknown nucleus '{}'. Supported: 1H, 13C, 19F, 31P, 15N, 11B, 29Si, 77Se, 17O, 33S",
1360 nucleus
1361 ))
1362 }
1363 };
1364 Ok(nmr::compute_nmr_spectrum(
1365 &shifts, &couplings, nuc, gamma, ppm_min, ppm_max, n_points,
1366 ))
1367}
1368
1369pub fn compute_hose_codes(smiles: &str, max_radius: usize) -> Result<Vec<nmr::HoseCode>, String> {
1371 let mol = graph::Molecule::from_smiles(smiles)?;
1372 Ok(nmr::hose::generate_hose_codes(&mol, max_radius))
1373}
1374
1375pub fn compute_orbital_mesh(
1387 elements: &[u8],
1388 positions: &[[f64; 3]],
1389 method: &str,
1390 mo_index: usize,
1391 spacing: f64,
1392 padding: f64,
1393 isovalue: f32,
1394) -> Result<mesh::OrbitalMeshResult, String> {
1395 let m = match method.to_lowercase().as_str() {
1396 "eht" | "huckel" => mesh::MeshMethod::Eht,
1397 "pm3" => mesh::MeshMethod::Pm3,
1398 "xtb" | "gfn0" | "gfn-xtb" => mesh::MeshMethod::Xtb,
1399 "hf3c" | "hf-3c" | "hf" => mesh::MeshMethod::Hf3c,
1400 _ => {
1401 return Err(format!(
1402 "Unknown method '{}'. Supported: eht, pm3, xtb, hf3c",
1403 method
1404 ))
1405 }
1406 };
1407 mesh::compute_orbital_mesh(elements, positions, m, mo_index, spacing, padding, isovalue)
1408}
1409
1410pub fn compute_ml_descriptors(
1417 elements: &[u8],
1418 bonds: &[(usize, usize, u8)],
1419 charges: &[f64],
1420 aromatic_atoms: &[bool],
1421) -> ml::MolecularDescriptors {
1422 ml::compute_descriptors(elements, bonds, charges, aromatic_atoms)
1423}
1424
1425pub fn predict_ml_properties(desc: &ml::MolecularDescriptors) -> ml::MlPropertyResult {
1430 ml::predict_properties(desc)
1431}
1432
1433pub fn compute_md_trajectory(
1435 smiles: &str,
1436 coords: &[f64],
1437 n_steps: usize,
1438 dt_fs: f64,
1439 seed: u64,
1440) -> Result<dynamics::MdTrajectory, String> {
1441 dynamics::simulate_velocity_verlet_uff(smiles, coords, n_steps, dt_fs, seed, None)
1442}
1443
1444pub fn compute_md_trajectory_nvt(
1446 smiles: &str,
1447 coords: &[f64],
1448 n_steps: usize,
1449 dt_fs: f64,
1450 seed: u64,
1451 target_temp_k: f64,
1452 thermostat_tau_fs: f64,
1453) -> Result<dynamics::MdTrajectory, String> {
1454 dynamics::simulate_velocity_verlet_uff(
1455 smiles,
1456 coords,
1457 n_steps,
1458 dt_fs,
1459 seed,
1460 Some((target_temp_k, thermostat_tau_fs)),
1461 )
1462}
1463
1464pub fn compute_simplified_neb_path(
1466 smiles: &str,
1467 start_coords: &[f64],
1468 end_coords: &[f64],
1469 n_images: usize,
1470 n_iter: usize,
1471 spring_k: f64,
1472 step_size: f64,
1473) -> Result<dynamics::NebPathResult, String> {
1474 dynamics::compute_simplified_neb_path(
1475 smiles,
1476 start_coords,
1477 end_coords,
1478 n_images,
1479 n_iter,
1480 spring_k,
1481 step_size,
1482 )
1483}
1484
1485fn coords_flat_to_matrix_f32(coords: &[f64], n_atoms: usize) -> nalgebra::DMatrix<f32> {
1486 let mut m = nalgebra::DMatrix::<f32>::zeros(n_atoms, 3);
1487 for i in 0..n_atoms {
1488 m[(i, 0)] = coords[3 * i] as f32;
1489 m[(i, 1)] = coords[3 * i + 1] as f32;
1490 m[(i, 2)] = coords[3 * i + 2] as f32;
1491 }
1492 m
1493}
1494
1495fn coords_matrix_f32_to_flat(m: &nalgebra::DMatrix<f32>) -> Vec<f64> {
1496 let mut out = Vec::with_capacity(m.nrows() * 3);
1497 for i in 0..m.nrows() {
1498 out.push(m[(i, 0)] as f64);
1499 out.push(m[(i, 1)] as f64);
1500 out.push(m[(i, 2)] as f64);
1501 }
1502 out
1503}
1504
1505pub fn search_conformers_with_uff(
1514 smiles: &str,
1515 n_samples: usize,
1516 seed: u64,
1517 rmsd_threshold: f64,
1518) -> Result<ConformerSearchResult, String> {
1519 if n_samples == 0 {
1520 return Err("n_samples must be > 0".to_string());
1521 }
1522 if rmsd_threshold <= 0.0 {
1523 return Err("rmsd_threshold must be > 0".to_string());
1524 }
1525
1526 let mol = graph::Molecule::from_smiles(smiles)?;
1527 let n_atoms = mol.graph.node_count();
1528 let bounds = distgeom::smooth_bounds_matrix(distgeom::calculate_bounds_matrix(&mol));
1529
1530 let mut generated = Vec::new();
1531 let mut notes = Vec::new();
1532 let mut rotatable_bonds = 0usize;
1533
1534 for i in 0..n_samples {
1535 let sample_seed = seed.wrapping_add(i as u64);
1536 let conf = embed(smiles, sample_seed);
1537
1538 if conf.error.is_some() || conf.coords.len() != n_atoms * 3 {
1539 continue;
1540 }
1541
1542 let mut coords = coords_flat_to_matrix_f32(&conf.coords, n_atoms);
1543 let rot_mc = forcefield::optimize_torsions_monte_carlo_bounds(
1544 &mut coords,
1545 &mol,
1546 &bounds,
1547 sample_seed ^ 0x9E37_79B9_7F4A_7C15,
1548 64,
1549 0.4,
1550 );
1551 let rot_greedy = forcefield::optimize_torsions_bounds(&mut coords, &mol, &bounds, 2);
1552 let rot = rot_mc.max(rot_greedy);
1553 rotatable_bonds = rot;
1554 let coords_flat = coords_matrix_f32_to_flat(&coords);
1555
1556 match compute_uff_energy(smiles, &coords_flat) {
1557 Ok(energy_kcal_mol) => generated.push(ConformerEnsembleMember {
1558 seed: sample_seed,
1559 cluster_id: None,
1560 coords: coords_flat,
1561 energy_kcal_mol,
1562 }),
1563 Err(_) => continue,
1564 }
1565 }
1566
1567 if generated.is_empty() {
1568 return Err("failed to generate any valid conformers".to_string());
1569 }
1570
1571 generated.sort_by(|a, b| {
1572 a.energy_kcal_mol
1573 .partial_cmp(&b.energy_kcal_mol)
1574 .unwrap_or(std::cmp::Ordering::Equal)
1575 });
1576
1577 let generated_count = generated.len();
1578
1579 let mut unique = Vec::new();
1580 let mut cluster_members: Vec<Vec<u64>> = Vec::new();
1581 for candidate in generated {
1582 let existing_cluster = unique.iter().position(|u: &ConformerEnsembleMember| {
1583 compute_rmsd(&candidate.coords, &u.coords) < rmsd_threshold
1584 });
1585
1586 if let Some(cluster_id) = existing_cluster {
1587 cluster_members[cluster_id].push(candidate.seed);
1588 } else {
1589 unique.push(candidate.clone());
1590 cluster_members.push(vec![candidate.seed]);
1591 }
1592 }
1593
1594 for (cluster_id, member) in unique.iter_mut().enumerate() {
1595 member.cluster_id = Some(cluster_id);
1596 }
1597
1598 let clusters: Vec<ConformerClusterSummary> = unique
1599 .iter()
1600 .enumerate()
1601 .map(|(cluster_id, representative)| ConformerClusterSummary {
1602 cluster_id,
1603 representative_seed: representative.seed,
1604 size: cluster_members[cluster_id].len(),
1605 member_seeds: cluster_members[cluster_id].clone(),
1606 })
1607 .collect();
1608
1609 notes.push(
1610 "Conformers are preconditioned with Monte Carlo torsion sampling + greedy torsion refinement, ranked by UFF energy, deduplicated by Kabsch-aligned RMSD threshold, and summarized as explicit RMSD clusters."
1611 .to_string(),
1612 );
1613
1614 Ok(ConformerSearchResult {
1615 generated: generated_count,
1616 unique: unique.len(),
1617 rotatable_bonds,
1618 conformers: unique,
1619 clusters,
1620 notes,
1621 })
1622}
1623
1624pub fn compute_uff_energy_with_aromatic_heuristics(
1626 smiles: &str,
1627 coords: &[f64],
1628) -> Result<reactivity::UffHeuristicEnergy, String> {
1629 let mol = graph::Molecule::from_smiles(smiles)?;
1630 let n = mol.graph.node_count();
1631 if coords.len() != n * 3 {
1632 return Err(format!("coords length {} != 3 * atoms {}", coords.len(), n));
1633 }
1634
1635 let ff = forcefield::builder::build_uff_force_field(&mol);
1636 let mut gradient = vec![0.0f64; n * 3];
1637 let raw = ff.compute_system_energy_and_gradients(coords, &mut gradient);
1638 Ok(reactivity::apply_aromatic_uff_correction(&mol, raw))
1639}
1640
1641pub fn compute_empirical_pka(smiles: &str) -> Result<reactivity::EmpiricalPkaResult, String> {
1643 let mol = graph::Molecule::from_smiles(smiles)?;
1644 let charges = compute_charges(smiles)?;
1645 Ok(reactivity::estimate_empirical_pka(&mol, &charges.charges))
1646}
1647
1648pub fn create_unit_cell(
1650 a: f64,
1651 b: f64,
1652 c: f64,
1653 alpha: f64,
1654 beta: f64,
1655 gamma: f64,
1656) -> materials::UnitCell {
1657 materials::UnitCell::from_parameters(&materials::CellParameters {
1658 a,
1659 b,
1660 c,
1661 alpha,
1662 beta,
1663 gamma,
1664 })
1665}
1666
1667pub fn assemble_framework(
1671 node: &materials::Sbu,
1672 linker: &materials::Sbu,
1673 topology: &materials::Topology,
1674 cell: &materials::UnitCell,
1675) -> materials::CrystalStructure {
1676 materials::assemble_framework(node, linker, topology, cell)
1677}
1678
1679pub fn compute_hf3c(
1687 elements: &[u8],
1688 positions: &[[f64; 3]],
1689 config: &hf::HfConfig,
1690) -> Result<hf::Hf3cResult, String> {
1691 hf::solve_hf3c(elements, positions, config)
1692}
1693
1694pub fn compute_ani(elements: &[u8], positions: &[[f64; 3]]) -> Result<ani::AniResult, String> {
1702 ani::api::compute_ani_test(elements, positions)
1703}
1704
1705pub fn compute_ani_with_models(
1712 elements: &[u8],
1713 positions: &[[f64; 3]],
1714 config: &ani::AniConfig,
1715 models: &std::collections::HashMap<u8, ani::nn::FeedForwardNet>,
1716) -> Result<ani::AniResult, String> {
1717 ani::compute_ani(elements, positions, config, models)
1718}
1719
1720pub fn compute_esp_grid(
1724 elements: &[u8],
1725 positions: &[[f64; 3]],
1726 spacing: f64,
1727 padding: f64,
1728) -> Result<esp::EspGrid, String> {
1729 compute_esp(elements, positions, spacing, padding)
1730}
1731
1732pub fn analyze_stereo(smiles: &str, coords: &[f64]) -> Result<stereo::StereoAnalysis, String> {
1736 let mol = graph::Molecule::from_smiles(smiles)?;
1737 let positions: Vec<[f64; 3]> = coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
1738 Ok(stereo::analyze_stereo(&mol, &positions))
1739}
1740
1741pub fn compute_nonpolar_solvation(
1745 elements: &[u8],
1746 positions: &[[f64; 3]],
1747 probe_radius: Option<f64>,
1748) -> solvation::NonPolarSolvation {
1749 solvation::compute_nonpolar_solvation(elements, positions, probe_radius)
1750}
1751
1752pub fn compute_gb_solvation(
1754 elements: &[u8],
1755 positions: &[[f64; 3]],
1756 charges: &[f64],
1757 solvent_dielectric: Option<f64>,
1758 solute_dielectric: Option<f64>,
1759 probe_radius: Option<f64>,
1760) -> solvation::GbSolvation {
1761 solvation::compute_gb_solvation(
1762 elements,
1763 positions,
1764 charges,
1765 solvent_dielectric,
1766 solute_dielectric,
1767 probe_radius,
1768 )
1769}
1770
1771pub fn butina_cluster(conformers: &[Vec<f64>], rmsd_cutoff: f64) -> clustering::ClusterResult {
1775 clustering::butina_cluster(conformers, rmsd_cutoff)
1776}
1777
1778pub fn compute_rmsd_matrix(conformers: &[Vec<f64>]) -> Vec<Vec<f64>> {
1780 clustering::compute_rmsd_matrix(conformers)
1781}
1782
1783pub fn compute_sssr(smiles: &str) -> Result<rings::sssr::SssrResult, String> {
1787 let mol = graph::Molecule::from_smiles(smiles)?;
1788 Ok(rings::sssr::compute_sssr(&mol))
1789}
1790
1791pub fn compute_ecfp(
1793 smiles: &str,
1794 radius: usize,
1795 n_bits: usize,
1796) -> Result<rings::ecfp::ECFPFingerprint, String> {
1797 let mol = graph::Molecule::from_smiles(smiles)?;
1798 Ok(rings::ecfp::compute_ecfp(&mol, radius, n_bits))
1799}
1800
1801pub fn compute_tanimoto(
1803 fp1: &rings::ecfp::ECFPFingerprint,
1804 fp2: &rings::ecfp::ECFPFingerprint,
1805) -> f64 {
1806 rings::ecfp::compute_tanimoto(fp1, fp2)
1807}
1808
1809pub fn compute_charges_configured(
1813 smiles: &str,
1814 config: &charges::gasteiger::GasteigerConfig,
1815) -> Result<charges::gasteiger::ChargeResult, String> {
1816 let mol = graph::Molecule::from_smiles(smiles)?;
1817 let n = mol.graph.node_count();
1818 let elements: Vec<u8> = (0..n)
1819 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
1820 .collect();
1821 let bonds: Vec<(usize, usize)> = mol
1822 .graph
1823 .edge_indices()
1824 .map(|e| {
1825 let (a, b) = mol.graph.edge_endpoints(e).unwrap();
1826 (a.index(), b.index())
1827 })
1828 .collect();
1829 let formal_charges: Vec<i8> = (0..n)
1830 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].formal_charge)
1831 .collect();
1832 charges::gasteiger::gasteiger_marsili_charges_configured(
1833 &elements,
1834 &bonds,
1835 &formal_charges,
1836 config,
1837 )
1838}
1839
1840pub fn compute_ensemble_j_couplings(
1844 smiles: &str,
1845 conformer_coords: &[Vec<f64>],
1846 energies_kcal: &[f64],
1847 temperature_k: f64,
1848) -> Result<Vec<nmr::JCoupling>, String> {
1849 let mol = graph::Molecule::from_smiles(smiles)?;
1850 let positions: Vec<Vec<[f64; 3]>> = conformer_coords
1851 .iter()
1852 .map(|c| c.chunks(3).map(|p| [p[0], p[1], p[2]]).collect())
1853 .collect();
1854 Ok(nmr::coupling::ensemble_averaged_j_couplings(
1855 &mol,
1856 &positions,
1857 energies_kcal,
1858 temperature_k,
1859 ))
1860}
1861
1862pub fn compute_ir_spectrum_broadened(
1866 analysis: &ir::VibrationalAnalysis,
1867 gamma: f64,
1868 wn_min: f64,
1869 wn_max: f64,
1870 n_points: usize,
1871 broadening: &str,
1872) -> ir::IrSpectrum {
1873 let bt = match broadening {
1874 "gaussian" | "Gaussian" => ir::BroadeningType::Gaussian,
1875 _ => ir::BroadeningType::Lorentzian,
1876 };
1877 ir::vibrations::compute_ir_spectrum_with_broadening(
1878 analysis, gamma, wn_min, wn_max, n_points, bt,
1879 )
1880}
1881
1882#[cfg(test)]
1883mod tests {
1884 use super::*;
1885
1886 #[test]
1887 fn test_cisplatin_style_support_metadata() {
1888 let smiles = "[Pt](Cl)(Cl)([NH3])[NH3]";
1889 let mol = parse(smiles).expect("Cisplatin-style example should parse");
1890 let elements: Vec<u8> = (0..mol.graph.node_count())
1891 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
1892 .collect();
1893
1894 let caps = get_system_capabilities(&elements);
1895 assert!(caps.eht.available);
1896 assert!(matches!(
1897 caps.eht.confidence,
1898 eht::SupportLevel::Experimental
1899 ));
1900 }
1901
1902 #[test]
1903 fn test_pt_system_routes_to_uff_when_experimental_disabled() {
1904 let smiles = "[Pt](Cl)(Cl)([NH3])[NH3]";
1905 let mol = parse(smiles).expect("Pt example should parse");
1906 let n = mol.graph.node_count();
1907 let elements: Vec<u8> = (0..n)
1908 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
1909 .collect();
1910 let positions = vec![[0.0, 0.0, 0.0]; n];
1911
1912 let result = compute_eht_or_uff_fallback(smiles, &elements, &positions, false)
1913 .expect("Fallback workflow should return a result");
1914
1915 assert!(matches!(
1916 result,
1917 ElectronicWorkflowResult::UffFallback { .. }
1918 ));
1919 }
1920
1921 #[test]
1922 fn test_method_plan_prefers_supported_workflows_for_organic_systems() {
1923 let plan = get_system_method_plan(&[6, 1, 1, 1, 1]);
1924
1925 assert_eq!(plan.geometry.recommended, Some(ScientificMethod::Embed));
1926 assert_eq!(
1927 plan.force_field_energy.recommended,
1928 Some(ScientificMethod::Uff)
1929 );
1930 assert_eq!(plan.orbitals.recommended, Some(ScientificMethod::Eht));
1931 assert_eq!(plan.population.recommended, Some(ScientificMethod::Eht));
1932 assert_eq!(plan.orbital_grid.recommended, Some(ScientificMethod::Eht));
1933 assert!(plan.orbitals.methods[0].confidence_score > 0.9);
1934 }
1935
1936 #[test]
1937 fn test_method_plan_marks_metal_orbital_workflow_experimental() {
1938 let plan = get_system_method_plan(&[78, 17, 17, 7, 7]);
1939
1940 assert_eq!(plan.orbitals.recommended, Some(ScientificMethod::Eht));
1941 assert_eq!(
1942 plan.force_field_energy.recommended,
1943 Some(ScientificMethod::Uff)
1944 );
1945 assert!(matches!(
1946 plan.orbitals.methods[0].confidence,
1947 eht::SupportLevel::Experimental
1948 ));
1949 assert!(plan.orbitals.methods[0].confidence_score < 0.9);
1950 assert!(!plan.orbitals.methods[0].warnings.is_empty());
1951 }
1952
1953 #[test]
1954 fn test_method_plan_reports_unavailable_workflows_for_unsupported_elements() {
1955 let plan = get_system_method_plan(&[92]);
1956
1957 assert_eq!(plan.force_field_energy.recommended, None);
1958 assert_eq!(plan.orbitals.recommended, None);
1959 assert_eq!(plan.population.recommended, None);
1960 assert_eq!(plan.orbital_grid.recommended, None);
1961 assert!(!plan.orbitals.methods[0].limitations.is_empty());
1962 }
1963
1964 #[test]
1965 fn test_compare_methods_supported_system_returns_success_rows() {
1966 let result = compare_methods("CC", &[6, 6], &[[0.0, 0.0, 0.0], [1.54, 0.0, 0.0]], true);
1967 assert_eq!(result.comparisons.len(), 2);
1968 assert!(result
1969 .comparisons
1970 .iter()
1971 .any(|entry| matches!(entry.method, ScientificMethod::Uff) && entry.available));
1972 assert!(result.comparisons.iter().any(|entry| matches!(
1973 entry.method,
1974 ScientificMethod::Eht
1975 ) && matches!(
1976 entry.status,
1977 MethodComparisonStatus::Success
1978 )));
1979 }
1980
1981 #[test]
1982 fn test_compare_methods_blocks_experimental_eht_when_disabled() {
1983 let result = compare_methods("[O]", &[78], &[[0.0, 0.0, 0.0]], false);
1984 let eht_row = result
1985 .comparisons
1986 .iter()
1987 .find(|entry| matches!(entry.method, ScientificMethod::Eht))
1988 .expect("EHT row must exist");
1989 assert!(matches!(
1990 eht_row.status,
1991 MethodComparisonStatus::Unavailable
1992 ));
1993 }
1994
1995 #[test]
1996 fn test_compare_methods_reports_unavailable_for_unsupported_elements() {
1997 let result = compare_methods("[O]", &[92], &[[0.0, 0.0, 0.0]], true);
1998 assert!(result
1999 .comparisons
2000 .iter()
2001 .all(|entry| matches!(entry.status, MethodComparisonStatus::Unavailable)));
2002 }
2003
2004 #[test]
2005 fn test_compute_fukui_descriptors_returns_atomwise_output() {
2006 let elements = [8u8, 1, 1];
2007 let positions = [[0.0, 0.0, 0.0], [0.757, 0.586, 0.0], [-0.757, 0.586, 0.0]];
2008 let result = compute_fukui_descriptors(&elements, &positions).unwrap();
2009 assert_eq!(result.num_atoms, 3);
2010 assert_eq!(result.condensed.len(), 3);
2011 }
2012
2013 #[test]
2014 fn test_compute_uv_vis_spectrum_returns_requested_grid_size() {
2015 let elements = [6u8, 6, 1, 1, 1, 1];
2016 let positions = [
2017 [0.0, 0.0, 0.0],
2018 [1.34, 0.0, 0.0],
2019 [-0.6, 0.92, 0.0],
2020 [-0.6, -0.92, 0.0],
2021 [1.94, 0.92, 0.0],
2022 [1.94, -0.92, 0.0],
2023 ];
2024 let spectrum = compute_uv_vis_spectrum(&elements, &positions, 0.2, 0.5, 8.0, 256).unwrap();
2025 assert_eq!(spectrum.energies_ev.len(), 256);
2026 assert_eq!(spectrum.intensities.len(), 256);
2027 }
2028
2029 #[test]
2030 fn test_analyze_graph_features_reports_benzene_aromaticity() {
2031 let analysis = analyze_graph_features("c1ccccc1").unwrap();
2032 assert_eq!(analysis.aromaticity.aromatic_atoms.len(), 12);
2033 assert_eq!(
2034 analysis
2035 .aromaticity
2036 .aromatic_atoms
2037 .iter()
2038 .filter(|v| **v)
2039 .count(),
2040 6
2041 );
2042 assert_eq!(analysis.aromaticity.aromatic_bonds.len(), 6);
2043 }
2044
2045 #[test]
2046 fn test_compute_empirical_pka_finds_acidic_site_for_acetic_acid() {
2047 let result = compute_empirical_pka("CC(=O)O").unwrap();
2048 assert!(!result.acidic_sites.is_empty());
2049 }
2050
2051 #[test]
2052 fn test_compute_uff_energy_with_aromatic_heuristics_applies_correction() {
2053 let conf = embed("c1ccccc1", 42);
2054 assert!(conf.error.is_none());
2055
2056 let result = compute_uff_energy_with_aromatic_heuristics("c1ccccc1", &conf.coords).unwrap();
2057 assert!(result.aromatic_bond_count >= 6);
2058 assert!(result.corrected_energy_kcal_mol <= result.raw_energy_kcal_mol);
2059 }
2060
2061 #[test]
2062 fn test_search_conformers_with_uff_returns_ranked_unique_ensemble() {
2063 let result = search_conformers_with_uff("CCCC", 10, 42, 0.2).unwrap();
2064 assert!(result.generated >= 1);
2065 assert!(result.unique >= 1);
2066 assert_eq!(result.unique, result.conformers.len());
2067 assert_eq!(result.unique, result.clusters.len());
2068 assert!(result.rotatable_bonds >= 1);
2069
2070 let mut total_members = 0usize;
2071 for (i, cluster) in result.clusters.iter().enumerate() {
2072 assert_eq!(cluster.cluster_id, i);
2073 assert!(cluster.size >= 1);
2074 total_members += cluster.size;
2075 assert_eq!(result.conformers[i].cluster_id, Some(i));
2076 assert_eq!(result.conformers[i].seed, cluster.representative_seed);
2077 }
2078 assert_eq!(total_members, result.generated);
2079
2080 for i in 1..result.conformers.len() {
2081 assert!(
2082 result.conformers[i - 1].energy_kcal_mol <= result.conformers[i].energy_kcal_mol
2083 );
2084 }
2085 }
2086
2087 #[test]
2088 fn test_search_conformers_with_uff_large_rmsd_threshold_collapses_duplicates() {
2089 let result = search_conformers_with_uff("CCCC", 8, 123, 10.0).unwrap();
2090 assert_eq!(result.unique, 1);
2091 assert_eq!(result.clusters.len(), 1);
2092 assert_eq!(result.clusters[0].size, result.generated);
2093 }
2094
2095 #[test]
2096 fn test_compute_md_trajectory_velocity_verlet_runs() {
2097 let conf = embed("CC", 42);
2098 assert!(conf.error.is_none());
2099
2100 let trj = compute_md_trajectory("CC", &conf.coords, 10, 0.25, 7).unwrap();
2101 assert_eq!(trj.frames.len(), 11);
2102 assert!(trj
2103 .frames
2104 .iter()
2105 .all(|f| f.coords.iter().all(|v| v.is_finite())));
2106 }
2107
2108 #[test]
2109 fn test_compute_md_trajectory_nvt_runs() {
2110 let conf = embed("CCO", 42);
2111 assert!(conf.error.is_none());
2112
2113 let trj =
2114 compute_md_trajectory_nvt("CCO", &conf.coords, 12, 0.25, 17, 300.0, 10.0).unwrap();
2115 assert_eq!(trj.frames.len(), 13);
2116 assert!(trj.frames.iter().all(|f| f.temperature_k.is_finite()));
2117 }
2118
2119 #[test]
2120 fn test_compute_simplified_neb_path_runs() {
2121 let c1 = embed("CC", 42);
2122 let c2 = embed("CC", 43);
2123 assert!(c1.error.is_none());
2124 assert!(c2.error.is_none());
2125
2126 let path =
2127 compute_simplified_neb_path("CC", &c1.coords, &c2.coords, 6, 20, 0.01, 1e-5).unwrap();
2128 assert_eq!(path.images.len(), 6);
2129 assert!(path
2130 .images
2131 .iter()
2132 .all(|img| img.potential_energy_kcal_mol.is_finite()));
2133 }
2134
2135 #[test]
2136 fn test_compute_hf3c_water() {
2137 let conf = embed("O", 42);
2138 assert!(conf.error.is_none());
2139 let pos: Vec<[f64; 3]> = conf.coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
2140
2141 let result = compute_hf3c(&conf.elements, &pos, &hf::HfConfig::default());
2142 assert!(result.is_ok(), "HF-3c should succeed for water");
2143 let r = result.unwrap();
2144 assert!(r.energy.is_finite());
2145 assert!(!r.orbital_energies.is_empty());
2146 }
2147
2148 #[test]
2149 fn test_compute_ani_water() {
2150 let conf = embed("O", 42);
2151 assert!(conf.error.is_none());
2152 let pos: Vec<[f64; 3]> = conf.coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
2153
2154 let result = compute_ani(&conf.elements, &pos);
2155 assert!(result.is_ok(), "ANI should succeed for water");
2156 let r = result.unwrap();
2157 assert!(r.energy.is_finite());
2158 assert_eq!(r.forces.len(), 3); assert_eq!(r.species.len(), 3);
2160 }
2161
2162 #[test]
2163 fn test_compute_esp_grid_water() {
2164 let conf = embed("O", 42);
2165 assert!(conf.error.is_none());
2166 let pos: Vec<[f64; 3]> = conf.coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
2167
2168 let result = compute_esp_grid(&conf.elements, &pos, 0.5, 3.0);
2169 assert!(result.is_ok(), "ESP grid should succeed for water");
2170 let g = result.unwrap();
2171 assert!(!g.values.is_empty());
2172 assert!(g.spacing > 0.0);
2173 assert!(g.dims[0] > 0 && g.dims[1] > 0 && g.dims[2] > 0);
2174 }
2175}