1#![allow(clippy::too_many_arguments)]
32#![allow(clippy::needless_range_loop)]
33
34pub mod alignment;
35pub mod ani;
36pub mod canonical_smiles;
37pub mod charges;
38pub mod charges_eeq;
39pub mod clustering;
40pub mod conformer;
41pub mod dipole;
42pub mod dispersion;
43pub mod distgeom;
44pub mod dos;
45pub mod dynamics;
46pub mod eht;
47pub mod esp;
48pub mod etkdg;
49pub mod experimental_status;
50#[doc(hidden)]
51pub mod fixture_io;
52pub mod forcefield;
53pub mod gpu;
54pub mod graph;
55pub mod hf;
56pub mod ir;
57pub mod materials;
58pub mod mesh;
59pub mod ml;
60pub mod nmr;
61pub mod optimization;
62pub mod periodic;
63pub mod pm3;
64pub mod population;
65pub mod reactivity;
66pub mod rings;
67pub mod scf;
68pub mod smarts;
69pub mod smiles;
70pub mod smirks;
71pub mod solvation;
72pub mod solvation_alpb;
73pub mod spectroscopy;
74pub mod stereo;
75pub mod surface;
76pub mod topology;
77pub mod transport;
78pub mod xtb;
79
80#[cfg(any(feature = "alpha-cga", feature = "alpha-gsm", feature = "alpha-sdr"))]
81pub mod alpha;
82
83#[cfg(any(
84 feature = "beta-kpm",
85 feature = "beta-mbh",
86 feature = "beta-cpm",
87 feature = "beta-randnla",
88 feature = "beta-riemannian"
89))]
90pub mod beta;
91
92#[cfg(feature = "alpha-dynamics-live")]
95pub mod dynamics_live;
96
97#[cfg(feature = "alpha-dft")]
98pub mod dft;
99
100#[cfg(feature = "alpha-mlff")]
101pub mod mlff {
102 pub use crate::ml::inference::*;
104 pub use crate::ml::mlff::*;
105 pub use crate::ml::symmetry_functions::*;
106}
107
108use serde::{Deserialize, Serialize};
109
110#[derive(Debug, Clone, Serialize, Deserialize)]
114pub struct ConformerResult {
115 pub smiles: String,
117 pub num_atoms: usize,
119 pub coords: Vec<f64>,
122 pub elements: Vec<u8>,
124 pub bonds: Vec<(usize, usize, String)>,
126 pub error: Option<String>,
128 pub time_ms: f64,
130}
131
132#[derive(Debug, Clone, Serialize, Deserialize)]
134pub struct ConformerConfig {
135 pub seed: u64,
137 pub num_threads: usize,
139}
140
141#[derive(Debug, Clone, Serialize, Deserialize)]
143pub struct ConformerEnsembleMember {
144 pub seed: u64,
146 pub cluster_id: Option<usize>,
148 pub coords: Vec<f64>,
150 pub energy_kcal_mol: f64,
152}
153
154#[derive(Debug, Clone, Serialize, Deserialize)]
156pub struct ConformerClusterSummary {
157 pub cluster_id: usize,
159 pub representative_seed: u64,
161 pub size: usize,
163 pub member_seeds: Vec<u64>,
165}
166
167#[derive(Debug, Clone, Serialize, Deserialize)]
169pub struct ConformerSearchResult {
170 pub generated: usize,
172 pub unique: usize,
174 pub rotatable_bonds: usize,
176 pub conformers: Vec<ConformerEnsembleMember>,
178 pub clusters: Vec<ConformerClusterSummary>,
180 pub notes: Vec<String>,
182}
183
184#[derive(Debug, Clone, Serialize, Deserialize)]
186pub struct MethodCapability {
187 pub available: bool,
189 pub confidence: eht::SupportLevel,
191 pub unsupported_elements: Vec<u8>,
193 pub warnings: Vec<String>,
195}
196
197#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
199#[serde(rename_all = "snake_case")]
200pub enum ScientificMethod {
201 Embed,
202 Uff,
203 Eht,
204 Pm3,
205 Xtb,
206 Mmff94,
207 Ani,
208 Hf3c,
209}
210
211#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
213#[serde(rename_all = "snake_case")]
214pub enum PropertyRequest {
215 Geometry,
216 ForceFieldEnergy,
217 Orbitals,
218 Population,
219 OrbitalGrid,
220}
221
222#[derive(Debug, Clone, Serialize, Deserialize)]
224pub struct MethodMetadata {
225 pub method: ScientificMethod,
226 pub available: bool,
227 pub confidence: eht::SupportLevel,
228 pub confidence_score: f64,
229 pub limitations: Vec<String>,
230 pub warnings: Vec<String>,
231}
232
233#[derive(Debug, Clone, Serialize, Deserialize)]
235pub struct PropertyMethodPlan {
236 pub property: PropertyRequest,
237 pub recommended: Option<ScientificMethod>,
238 pub fallback: Option<ScientificMethod>,
239 pub rationale: Vec<String>,
240 pub methods: Vec<MethodMetadata>,
241}
242
243#[derive(Debug, Clone, Serialize, Deserialize)]
245pub struct SystemMethodPlan {
246 pub capabilities: SystemCapabilities,
247 pub geometry: PropertyMethodPlan,
248 pub force_field_energy: PropertyMethodPlan,
249 pub orbitals: PropertyMethodPlan,
250 pub population: PropertyMethodPlan,
251 pub orbital_grid: PropertyMethodPlan,
252}
253
254#[derive(Debug, Clone, Serialize, Deserialize)]
256pub struct SystemCapabilities {
257 pub embed: MethodCapability,
258 pub uff: MethodCapability,
259 pub eht: MethodCapability,
260 pub population: MethodCapability,
261 pub orbital_grid: MethodCapability,
262}
263
264#[derive(Debug, Clone, Serialize, Deserialize)]
266#[serde(tag = "mode", rename_all = "snake_case")]
267pub enum ElectronicWorkflowResult {
268 Eht {
269 result: eht::EhtResult,
270 },
271 UffFallback {
272 energy_kcal_mol: f64,
273 reason: String,
274 support: eht::EhtSupport,
275 },
276}
277
278#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
280#[serde(rename_all = "snake_case")]
281pub enum MethodComparisonStatus {
282 Success,
283 Unavailable,
284 Error,
285}
286
287#[derive(Debug, Clone, Serialize, Deserialize)]
289#[serde(tag = "kind", rename_all = "snake_case")]
290pub enum MethodComparisonPayload {
291 Eht {
292 homo_energy: f64,
293 lumo_energy: f64,
294 gap: f64,
295 support: eht::EhtSupport,
296 },
297 Uff {
298 energy_kcal_mol: f64,
299 },
300}
301
302#[derive(Debug, Clone, Serialize, Deserialize)]
304pub struct MethodComparisonEntry {
305 pub method: ScientificMethod,
306 pub status: MethodComparisonStatus,
307 pub available: bool,
308 pub confidence: eht::SupportLevel,
309 pub confidence_score: f64,
310 pub warnings: Vec<String>,
311 pub limitations: Vec<String>,
312 pub payload: Option<MethodComparisonPayload>,
313 pub error: Option<String>,
314}
315
316#[derive(Debug, Clone, Serialize, Deserialize)]
318pub struct MethodComparisonResult {
319 pub plan: SystemMethodPlan,
320 pub comparisons: Vec<MethodComparisonEntry>,
321}
322
323#[derive(Debug, Clone, Serialize, Deserialize)]
325pub struct AromaticityAnalysis {
326 pub aromatic_atoms: Vec<bool>,
327 pub aromatic_bonds: Vec<(usize, usize)>,
328}
329
330#[derive(Debug, Clone, Serialize, Deserialize)]
332pub struct StereocenterAnalysis {
333 pub tagged_tetrahedral_centers: Vec<usize>,
334 pub inferred_tetrahedral_centers: Vec<usize>,
335}
336
337#[derive(Debug, Clone, Serialize, Deserialize)]
339pub struct GraphFeatureAnalysis {
340 pub aromaticity: AromaticityAnalysis,
341 pub stereocenters: StereocenterAnalysis,
342}
343
344#[derive(Debug, Clone, Serialize, Deserialize)]
346pub struct GiaoNmrConfig {
347 pub charge: i32,
349 pub multiplicity: u32,
351 pub max_scf_iterations: usize,
353 pub use_parallel_eri: bool,
355 pub allow_basis_fallback: bool,
357}
358
359#[derive(Debug, Clone, Serialize, Deserialize)]
361pub struct GiaoNmrEntry {
362 pub atom_index: usize,
363 pub element: u8,
364 pub tensor: [[f64; 3]; 3],
365 pub isotropic: f64,
366 pub anisotropy: f64,
367 pub chemical_shift: f64,
368}
369
370#[derive(Debug, Clone, Serialize, Deserialize)]
372pub struct GiaoNmrResult {
373 pub nucleus: String,
374 pub target_atomic_number: u8,
375 pub method: String,
376 pub basis_set: String,
377 pub charge: i32,
378 pub multiplicity: u32,
379 pub scf_converged: bool,
380 pub scf_iterations: usize,
381 pub total_energy_hartree: f64,
382 pub n_basis: usize,
383 pub n_target_atoms: usize,
384 pub chemical_shifts: Vec<f64>,
385 pub shieldings: Vec<GiaoNmrEntry>,
386 pub fallback_elements: Vec<u8>,
387 pub notes: Vec<String>,
388}
389
390impl Default for ConformerConfig {
391 fn default() -> Self {
392 Self {
393 seed: 42,
394 num_threads: 0,
395 }
396 }
397}
398
399impl Default for GiaoNmrConfig {
400 fn default() -> Self {
401 Self {
402 charge: 0,
403 multiplicity: 1,
404 max_scf_iterations: 100,
405 use_parallel_eri: false,
406 allow_basis_fallback: false,
407 }
408 }
409}
410
411pub fn version() -> String {
415 format!("sci-form {}", env!("CARGO_PKG_VERSION"))
416}
417
418pub fn get_eht_support(elements: &[u8]) -> eht::EhtSupport {
420 eht::analyze_eht_support(elements)
421}
422
423fn is_uff_element_supported(z: u8) -> bool {
424 matches!(
425 z,
426 1 | 5
427 | 6
428 | 7
429 | 8
430 | 9
431 | 14
432 | 15
433 | 16
434 | 17
435 | 22
436 | 23
437 | 24
438 | 25
439 | 26
440 | 27
441 | 28
442 | 29
443 | 30
444 | 32
445 | 33
446 | 34
447 | 35
448 | 42
449 | 46
450 | 47
451 | 50
452 | 51
453 | 52
454 | 53
455 | 78
456 | 79
457 )
458}
459
460fn unique_sorted_unsupported(elements: &[u8], pred: impl Fn(u8) -> bool) -> Vec<u8> {
461 let mut out: Vec<u8> = elements.iter().copied().filter(|&z| !pred(z)).collect();
462 out.sort_unstable();
463 out.dedup();
464 out
465}
466
467pub fn get_system_capabilities(elements: &[u8]) -> SystemCapabilities {
469 let eht_support = get_eht_support(elements);
470 let uff_unsupported = unique_sorted_unsupported(elements, is_uff_element_supported);
471
472 let embed = MethodCapability {
473 available: !elements.is_empty(),
474 confidence: eht::SupportLevel::Experimental,
475 unsupported_elements: Vec::new(),
476 warnings: vec![
477 "Embed capability is inferred from element presence only; final success still depends on full molecular graph and geometry constraints.".to_string(),
478 ],
479 };
480
481 let uff = if uff_unsupported.is_empty() {
482 MethodCapability {
483 available: true,
484 confidence: eht::SupportLevel::Supported,
485 unsupported_elements: Vec::new(),
486 warnings: Vec::new(),
487 }
488 } else {
489 MethodCapability {
490 available: false,
491 confidence: eht::SupportLevel::Unsupported,
492 unsupported_elements: uff_unsupported.clone(),
493 warnings: vec![format!(
494 "UFF atom typing is unavailable for elements {:?}.",
495 uff_unsupported
496 )],
497 }
498 };
499
500 let eht = MethodCapability {
501 available: eht_support.unsupported_elements.is_empty(),
502 confidence: eht_support.level,
503 unsupported_elements: eht_support.unsupported_elements.clone(),
504 warnings: eht_support.warnings.clone(),
505 };
506
507 let population = MethodCapability {
508 available: eht.available,
509 confidence: eht.confidence,
510 unsupported_elements: eht.unsupported_elements.clone(),
511 warnings: eht.warnings.clone(),
512 };
513
514 let orbital_grid = MethodCapability {
515 available: eht.available,
516 confidence: eht.confidence,
517 unsupported_elements: eht.unsupported_elements.clone(),
518 warnings: eht.warnings.clone(),
519 };
520
521 SystemCapabilities {
522 embed,
523 uff,
524 eht,
525 population,
526 orbital_grid,
527 }
528}
529
530fn confidence_score_for_method(method: ScientificMethod, capability: &MethodCapability) -> f64 {
531 if !capability.available {
532 return 0.0;
533 }
534
535 match method {
536 ScientificMethod::Embed => 0.8,
537 ScientificMethod::Uff | ScientificMethod::Mmff94 => match capability.confidence {
538 eht::SupportLevel::Supported => 0.95,
539 eht::SupportLevel::Experimental => 0.75,
540 eht::SupportLevel::Unsupported => 0.0,
541 },
542 ScientificMethod::Eht | ScientificMethod::Pm3 | ScientificMethod::Xtb => {
543 match capability.confidence {
544 eht::SupportLevel::Supported => 0.95,
545 eht::SupportLevel::Experimental => 0.6,
546 eht::SupportLevel::Unsupported => 0.0,
547 }
548 }
549 ScientificMethod::Ani => match capability.confidence {
550 eht::SupportLevel::Supported => 0.90,
551 eht::SupportLevel::Experimental => 0.7,
552 eht::SupportLevel::Unsupported => 0.0,
553 },
554 ScientificMethod::Hf3c => match capability.confidence {
555 eht::SupportLevel::Supported => 0.85,
556 eht::SupportLevel::Experimental => 0.65,
557 eht::SupportLevel::Unsupported => 0.0,
558 },
559 }
560}
561
562fn build_method_metadata(
563 method: ScientificMethod,
564 capability: &MethodCapability,
565 extra_limitations: &[&str],
566) -> MethodMetadata {
567 let mut limitations: Vec<String> = extra_limitations.iter().map(|s| s.to_string()).collect();
568
569 if !capability.unsupported_elements.is_empty() {
570 limitations.push(format!(
571 "Unsupported elements for this method: {:?}.",
572 capability.unsupported_elements
573 ));
574 }
575
576 if matches!(method, ScientificMethod::Eht)
577 && matches!(capability.confidence, eht::SupportLevel::Experimental)
578 {
579 limitations.push(
580 "Transition-metal EHT parameters remain provisional and should be treated as experimental."
581 .to_string(),
582 );
583 }
584
585 MethodMetadata {
586 method,
587 available: capability.available,
588 confidence: capability.confidence,
589 confidence_score: confidence_score_for_method(method, capability),
590 limitations,
591 warnings: capability.warnings.clone(),
592 }
593}
594
595fn build_property_plan(
596 property: PropertyRequest,
597 recommended: Option<ScientificMethod>,
598 fallback: Option<ScientificMethod>,
599 rationale: Vec<String>,
600 methods: Vec<MethodMetadata>,
601) -> PropertyMethodPlan {
602 PropertyMethodPlan {
603 property,
604 recommended,
605 fallback,
606 rationale,
607 methods,
608 }
609}
610
611pub fn get_system_method_plan(elements: &[u8]) -> SystemMethodPlan {
613 let capabilities = get_system_capabilities(elements);
614
615 let geometry_method = build_method_metadata(
616 ScientificMethod::Embed,
617 &capabilities.embed,
618 &["Geometry generation still depends on graph topology, stereochemistry, and embedding constraints."],
619 );
620 let geometry = build_property_plan(
621 PropertyRequest::Geometry,
622 geometry_method.available.then_some(ScientificMethod::Embed),
623 None,
624 vec!["Embedding is the top-level geometry generation path in sci-form.".to_string()],
625 vec![geometry_method],
626 );
627
628 let uff_method = build_method_metadata(
629 ScientificMethod::Uff,
630 &capabilities.uff,
631 &["This recommendation applies to force-field energy evaluation, not molecular orbital analysis."],
632 );
633 let force_field_energy = build_property_plan(
634 PropertyRequest::ForceFieldEnergy,
635 uff_method.available.then_some(ScientificMethod::Uff),
636 None,
637 vec![
638 "UFF is the top-level force-field energy path when atom typing is available."
639 .to_string(),
640 ],
641 vec![uff_method],
642 );
643
644 let eht_method = build_method_metadata(
645 ScientificMethod::Eht,
646 &capabilities.eht,
647 &["EHT is the only current top-level orbital method in sci-form."],
648 );
649 let orbitals = build_property_plan(
650 PropertyRequest::Orbitals,
651 eht_method.available.then_some(ScientificMethod::Eht),
652 None,
653 vec!["Orbital energies and MO coefficients are produced by the EHT workflow.".to_string()],
654 vec![eht_method.clone()],
655 );
656
657 let population_method = build_method_metadata(
658 ScientificMethod::Eht,
659 &capabilities.population,
660 &["Population analysis is derived from the EHT density and overlap matrices."],
661 );
662 let population = build_property_plan(
663 PropertyRequest::Population,
664 population_method.available.then_some(ScientificMethod::Eht),
665 None,
666 vec!["Population analysis currently requires a successful EHT calculation.".to_string()],
667 vec![population_method],
668 );
669
670 let orbital_grid_method = build_method_metadata(
671 ScientificMethod::Eht,
672 &capabilities.orbital_grid,
673 &["Orbital-grid rendering currently depends on EHT molecular-orbital coefficients."],
674 );
675 let orbital_grid = build_property_plan(
676 PropertyRequest::OrbitalGrid,
677 orbital_grid_method
678 .available
679 .then_some(ScientificMethod::Eht),
680 None,
681 vec![
682 "Orbital-grid generation currently requires a successful EHT calculation.".to_string(),
683 ],
684 vec![orbital_grid_method],
685 );
686
687 SystemMethodPlan {
688 capabilities,
689 geometry,
690 force_field_energy,
691 orbitals,
692 population,
693 orbital_grid,
694 }
695}
696
697pub fn embed(smiles: &str, seed: u64) -> ConformerResult {
699 #[cfg(not(target_arch = "wasm32"))]
700 let start = std::time::Instant::now();
701
702 let mol = match graph::Molecule::from_smiles(smiles) {
703 Ok(m) => m,
704 Err(e) => {
705 return ConformerResult {
706 smiles: smiles.to_string(),
707 num_atoms: 0,
708 coords: vec![],
709 elements: vec![],
710 bonds: vec![],
711 error: Some(e),
712 #[cfg(not(target_arch = "wasm32"))]
713 time_ms: start.elapsed().as_secs_f64() * 1000.0,
714 #[cfg(target_arch = "wasm32")]
715 time_ms: 0.0,
716 };
717 }
718 };
719
720 let n = mol.graph.node_count();
721 let elements: Vec<u8> = (0..n)
722 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
723 .collect();
724 let bonds: Vec<(usize, usize, String)> = mol
725 .graph
726 .edge_indices()
727 .map(|e| {
728 let (a, b) = mol.graph.edge_endpoints(e).unwrap();
729 let order = match mol.graph[e].order {
730 graph::BondOrder::Single => "SINGLE",
731 graph::BondOrder::Double => "DOUBLE",
732 graph::BondOrder::Triple => "TRIPLE",
733 graph::BondOrder::Aromatic => "AROMATIC",
734 graph::BondOrder::Unknown => "UNKNOWN",
735 };
736 (a.index(), b.index(), order.to_string())
737 })
738 .collect();
739
740 let max_seed_retries = 3u64;
743 let mut last_err = String::new();
744 for retry in 0..max_seed_retries {
745 let current_seed = seed.wrapping_add(retry.wrapping_mul(997));
746 match conformer::generate_3d_conformer(&mol, current_seed) {
747 Ok(coords) => {
748 let mut flat = Vec::with_capacity(n * 3);
749 for i in 0..n {
750 flat.push(coords[(i, 0)] as f64);
751 flat.push(coords[(i, 1)] as f64);
752 flat.push(coords[(i, 2)] as f64);
753 }
754 return ConformerResult {
755 smiles: smiles.to_string(),
756 num_atoms: n,
757 coords: flat,
758 elements,
759 bonds,
760 error: None,
761 #[cfg(not(target_arch = "wasm32"))]
762 time_ms: start.elapsed().as_secs_f64() * 1000.0,
763 #[cfg(target_arch = "wasm32")]
764 time_ms: 0.0,
765 };
766 }
767 Err(e) => {
768 last_err = e;
769 }
770 }
771 }
772 ConformerResult {
773 smiles: smiles.to_string(),
774 num_atoms: n,
775 coords: vec![],
776 elements,
777 bonds,
778 error: Some(last_err),
779 #[cfg(not(target_arch = "wasm32"))]
780 time_ms: start.elapsed().as_secs_f64() * 1000.0,
781 #[cfg(target_arch = "wasm32")]
782 time_ms: 0.0,
783 }
784}
785
786#[cfg(feature = "parallel")]
791pub fn embed_batch(smiles_list: &[&str], config: &ConformerConfig) -> Vec<ConformerResult> {
792 use rayon::prelude::*;
793
794 if config.num_threads > 0 {
795 let pool = rayon::ThreadPoolBuilder::new()
796 .num_threads(config.num_threads)
797 .build()
798 .unwrap();
799 pool.install(|| {
800 smiles_list
801 .par_iter()
802 .map(|smi| embed(smi, config.seed))
803 .collect()
804 })
805 } else {
806 smiles_list
807 .par_iter()
808 .map(|smi| embed(smi, config.seed))
809 .collect()
810 }
811}
812
813#[cfg(not(feature = "parallel"))]
815pub fn embed_batch(smiles_list: &[&str], config: &ConformerConfig) -> Vec<ConformerResult> {
816 smiles_list
817 .iter()
818 .map(|smi| embed(smi, config.seed))
819 .collect()
820}
821
822pub fn embed_diverse(
833 smiles: &str,
834 n_conformers: usize,
835 rmsd_cutoff: f64,
836 base_seed: u64,
837) -> Vec<ConformerResult> {
838 let all_results: Vec<ConformerResult> = (0..n_conformers as u64)
839 .map(|i| embed(smiles, base_seed.wrapping_add(i)))
840 .collect();
841
842 let successful: Vec<(usize, &ConformerResult)> = all_results
843 .iter()
844 .enumerate()
845 .filter(|(_, r)| r.error.is_none() && !r.coords.is_empty())
846 .collect();
847
848 if successful.len() <= 1 {
849 return all_results
850 .into_iter()
851 .filter(|r| r.error.is_none())
852 .collect();
853 }
854
855 let coords_vecs: Vec<Vec<f64>> = successful.iter().map(|(_, r)| r.coords.clone()).collect();
856 let cluster_result = clustering::butina_cluster(&coords_vecs, rmsd_cutoff);
857
858 cluster_result
859 .centroid_indices
860 .iter()
861 .map(|&ci| {
862 let (orig_idx, _) = successful[ci];
863 all_results[orig_idx].clone()
864 })
865 .collect()
866}
867
868pub fn parse(smiles: &str) -> Result<graph::Molecule, String> {
870 graph::Molecule::from_smiles(smiles)
871}
872
873pub fn compute_charges(smiles: &str) -> Result<charges::gasteiger::ChargeResult, String> {
878 let mol = graph::Molecule::from_smiles(smiles)?;
879 let n = mol.graph.node_count();
880 let elements: Vec<u8> = (0..n)
881 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
882 .collect();
883 let formal_charges: Vec<i8> = (0..n)
884 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].formal_charge)
885 .collect();
886 let bonds: Vec<(usize, usize)> = mol
887 .graph
888 .edge_indices()
889 .map(|e| {
890 let (a, b) = mol.graph.edge_endpoints(e).unwrap();
891 (a.index(), b.index())
892 })
893 .collect();
894 charges::gasteiger::gasteiger_marsili_charges(&elements, &bonds, &formal_charges, 6)
895}
896
897pub fn compute_sasa(
901 elements: &[u8],
902 coords_flat: &[f64],
903 probe_radius: Option<f64>,
904) -> Result<surface::sasa::SasaResult, String> {
905 if coords_flat.len() != elements.len() * 3 {
906 return Err(format!(
907 "coords length {} != 3 * elements {}",
908 coords_flat.len(),
909 elements.len()
910 ));
911 }
912 let positions: Vec<[f64; 3]> = coords_flat
913 .chunks_exact(3)
914 .map(|c| [c[0], c[1], c[2]])
915 .collect();
916
917 #[cfg(feature = "parallel")]
918 {
919 Ok(surface::sasa::compute_sasa_parallel(
920 elements,
921 &positions,
922 probe_radius,
923 None,
924 ))
925 }
926
927 #[cfg(not(feature = "parallel"))]
928 {
929 Ok(surface::sasa::compute_sasa(
930 elements,
931 &positions,
932 probe_radius,
933 None,
934 ))
935 }
936}
937
938pub fn compute_population(
940 elements: &[u8],
941 positions: &[[f64; 3]],
942) -> Result<population::PopulationResult, String> {
943 let eht_result = eht::solve_eht(elements, positions, None)?;
944 Ok(population::compute_population(
945 elements,
946 positions,
947 &eht_result.coefficients,
948 eht_result.n_electrons,
949 ))
950}
951
952pub fn compute_dipole(
954 elements: &[u8],
955 positions: &[[f64; 3]],
956) -> Result<dipole::DipoleResult, String> {
957 let eht_result = eht::solve_eht(elements, positions, None)?;
958 Ok(dipole::compute_dipole_from_eht(
959 elements,
960 positions,
961 &eht_result.coefficients,
962 eht_result.n_electrons,
963 ))
964}
965
966pub fn compute_frontier_descriptors(
968 elements: &[u8],
969 positions: &[[f64; 3]],
970) -> Result<reactivity::FrontierDescriptors, String> {
971 let eht_result = eht::solve_eht(elements, positions, None)?;
972 Ok(reactivity::compute_frontier_descriptors(
973 elements,
974 positions,
975 &eht_result,
976 ))
977}
978
979pub fn compute_fukui_descriptors(
981 elements: &[u8],
982 positions: &[[f64; 3]],
983) -> Result<reactivity::FukuiDescriptors, String> {
984 let eht_result = eht::solve_eht(elements, positions, None)?;
985 Ok(reactivity::compute_fukui_descriptors(
986 elements,
987 positions,
988 &eht_result,
989 ))
990}
991
992pub fn compute_reactivity_ranking(
994 elements: &[u8],
995 positions: &[[f64; 3]],
996) -> Result<reactivity::ReactivityRanking, String> {
997 let eht_result = eht::solve_eht(elements, positions, None)?;
998 let fukui = reactivity::compute_fukui_descriptors(elements, positions, &eht_result);
999 let pop = population::compute_population(
1000 elements,
1001 positions,
1002 &eht_result.coefficients,
1003 eht_result.n_electrons,
1004 );
1005 Ok(reactivity::rank_reactivity_sites(
1006 &fukui,
1007 &pop.mulliken_charges,
1008 ))
1009}
1010
1011pub fn compute_uv_vis_spectrum(
1013 elements: &[u8],
1014 positions: &[[f64; 3]],
1015 sigma: f64,
1016 e_min: f64,
1017 e_max: f64,
1018 n_points: usize,
1019) -> Result<reactivity::UvVisSpectrum, String> {
1020 let eht_result = eht::solve_eht(elements, positions, None)?;
1021 Ok(reactivity::compute_uv_vis_like_spectrum(
1022 &eht_result,
1023 sigma,
1024 e_min,
1025 e_max,
1026 n_points,
1027 ))
1028}
1029
1030pub fn compute_bond_orders(
1032 elements: &[u8],
1033 positions: &[[f64; 3]],
1034) -> Result<population::BondOrderResult, String> {
1035 let eht_result = eht::solve_eht(elements, positions, None)?;
1036 Ok(population::compute_bond_orders(
1037 elements,
1038 positions,
1039 &eht_result.coefficients,
1040 eht_result.n_electrons,
1041 ))
1042}
1043
1044pub fn compute_topology(
1046 elements: &[u8],
1047 positions: &[[f64; 3]],
1048) -> topology::TopologyAnalysisResult {
1049 topology::analyze_topology(elements, positions)
1050}
1051
1052pub fn analyze_graph_features(smiles: &str) -> Result<GraphFeatureAnalysis, String> {
1054 use petgraph::visit::EdgeRef;
1055
1056 let mol = parse(smiles)?;
1057 let n_atoms = mol.graph.node_count();
1058 let mut aromatic_atoms = vec![false; n_atoms];
1059 let mut aromatic_bonds = Vec::new();
1060
1061 for edge in mol.graph.edge_references() {
1062 if matches!(edge.weight().order, graph::BondOrder::Aromatic) {
1063 let i = edge.source().index();
1064 let j = edge.target().index();
1065 aromatic_atoms[i] = true;
1066 aromatic_atoms[j] = true;
1067 aromatic_bonds.push((i, j));
1068 }
1069 }
1070
1071 let mut tagged_tetrahedral_centers = Vec::new();
1072 let mut inferred_tetrahedral_centers = Vec::new();
1073 for i in 0..n_atoms {
1074 let idx = petgraph::graph::NodeIndex::new(i);
1075 let atom = &mol.graph[idx];
1076 if matches!(
1077 atom.chiral_tag,
1078 graph::ChiralType::TetrahedralCW | graph::ChiralType::TetrahedralCCW
1079 ) {
1080 tagged_tetrahedral_centers.push(i);
1081 }
1082
1083 let neighbors: Vec<_> = mol.graph.neighbors(idx).collect();
1084 if neighbors.len() == 4 && matches!(atom.hybridization, graph::Hybridization::SP3) {
1085 let mut signature: Vec<u8> = neighbors.iter().map(|n| mol.graph[*n].element).collect();
1086 signature.sort_unstable();
1087 signature.dedup();
1088 if signature.len() >= 3 {
1089 inferred_tetrahedral_centers.push(i);
1090 }
1091 }
1092 }
1093
1094 Ok(GraphFeatureAnalysis {
1095 aromaticity: AromaticityAnalysis {
1096 aromatic_atoms,
1097 aromatic_bonds,
1098 },
1099 stereocenters: StereocenterAnalysis {
1100 tagged_tetrahedral_centers,
1101 inferred_tetrahedral_centers,
1102 },
1103 })
1104}
1105
1106pub fn compute_eht_or_uff_fallback(
1111 smiles: &str,
1112 elements: &[u8],
1113 positions: &[[f64; 3]],
1114 allow_experimental_eht: bool,
1115) -> Result<ElectronicWorkflowResult, String> {
1116 let support = get_eht_support(elements);
1117 let should_fallback = match support.level {
1118 eht::SupportLevel::Unsupported => true,
1119 eht::SupportLevel::Experimental => !allow_experimental_eht,
1120 eht::SupportLevel::Supported => false,
1121 };
1122
1123 if should_fallback {
1124 let coords_flat: Vec<f64> = positions.iter().flat_map(|p| p.iter().copied()).collect();
1125 let energy = compute_uff_energy(smiles, &coords_flat).map_err(|e| {
1126 format!(
1127 "EHT is not appropriate for this system (support: {:?}) and UFF fallback failed: {}",
1128 support.level, e
1129 )
1130 })?;
1131 return Ok(ElectronicWorkflowResult::UffFallback {
1132 energy_kcal_mol: energy,
1133 reason: if matches!(support.level, eht::SupportLevel::Unsupported) {
1134 "EHT unsupported for one or more elements; routed to UFF-only workflow.".to_string()
1135 } else {
1136 "EHT confidence is experimental and experimental mode is disabled; routed to UFF-only workflow."
1137 .to_string()
1138 },
1139 support,
1140 });
1141 }
1142
1143 let result = eht::solve_eht(elements, positions, None)?;
1144 Ok(ElectronicWorkflowResult::Eht { result })
1145}
1146
1147pub fn compare_methods(
1152 smiles: &str,
1153 elements: &[u8],
1154 positions: &[[f64; 3]],
1155 allow_experimental_eht: bool,
1156) -> MethodComparisonResult {
1157 let plan = get_system_method_plan(elements);
1158 let mut comparisons = Vec::new();
1159
1160 let coords_flat: Vec<f64> = positions.iter().flat_map(|p| p.iter().copied()).collect();
1161
1162 {
1163 let meta = build_method_metadata(
1164 ScientificMethod::Uff,
1165 &plan.capabilities.uff,
1166 &["Comparison uses UFF force-field energy as the UFF observable."],
1167 );
1168 if !meta.available {
1169 comparisons.push(MethodComparisonEntry {
1170 method: ScientificMethod::Uff,
1171 status: MethodComparisonStatus::Unavailable,
1172 available: false,
1173 confidence: meta.confidence,
1174 confidence_score: meta.confidence_score,
1175 warnings: meta.warnings,
1176 limitations: meta.limitations,
1177 payload: None,
1178 error: Some("UFF is unavailable for this element set.".to_string()),
1179 });
1180 } else {
1181 match compute_uff_energy(smiles, &coords_flat) {
1182 Ok(energy) => comparisons.push(MethodComparisonEntry {
1183 method: ScientificMethod::Uff,
1184 status: MethodComparisonStatus::Success,
1185 available: true,
1186 confidence: meta.confidence,
1187 confidence_score: meta.confidence_score,
1188 warnings: meta.warnings,
1189 limitations: meta.limitations,
1190 payload: Some(MethodComparisonPayload::Uff {
1191 energy_kcal_mol: energy,
1192 }),
1193 error: None,
1194 }),
1195 Err(err) => comparisons.push(MethodComparisonEntry {
1196 method: ScientificMethod::Uff,
1197 status: MethodComparisonStatus::Error,
1198 available: true,
1199 confidence: meta.confidence,
1200 confidence_score: meta.confidence_score,
1201 warnings: meta.warnings,
1202 limitations: meta.limitations,
1203 payload: None,
1204 error: Some(err),
1205 }),
1206 }
1207 }
1208 }
1209
1210 {
1211 let meta = build_method_metadata(
1212 ScientificMethod::Eht,
1213 &plan.capabilities.eht,
1214 &["Comparison uses frontier orbital energies and gap as the EHT observable."],
1215 );
1216
1217 if !meta.available {
1218 comparisons.push(MethodComparisonEntry {
1219 method: ScientificMethod::Eht,
1220 status: MethodComparisonStatus::Unavailable,
1221 available: false,
1222 confidence: meta.confidence,
1223 confidence_score: meta.confidence_score,
1224 warnings: meta.warnings,
1225 limitations: meta.limitations,
1226 payload: None,
1227 error: Some("EHT is unavailable for this element set.".to_string()),
1228 });
1229 } else if matches!(meta.confidence, eht::SupportLevel::Experimental)
1230 && !allow_experimental_eht
1231 {
1232 comparisons.push(MethodComparisonEntry {
1233 method: ScientificMethod::Eht,
1234 status: MethodComparisonStatus::Unavailable,
1235 available: true,
1236 confidence: meta.confidence,
1237 confidence_score: meta.confidence_score,
1238 warnings: meta.warnings,
1239 limitations: meta.limitations,
1240 payload: None,
1241 error: Some(
1242 "EHT confidence is experimental and allow_experimental_eht=false.".to_string(),
1243 ),
1244 });
1245 } else {
1246 match eht::solve_eht(elements, positions, None) {
1247 Ok(result) => comparisons.push(MethodComparisonEntry {
1248 method: ScientificMethod::Eht,
1249 status: MethodComparisonStatus::Success,
1250 available: true,
1251 confidence: meta.confidence,
1252 confidence_score: meta.confidence_score,
1253 warnings: meta.warnings,
1254 limitations: meta.limitations,
1255 payload: Some(MethodComparisonPayload::Eht {
1256 homo_energy: result.homo_energy,
1257 lumo_energy: result.lumo_energy,
1258 gap: result.gap,
1259 support: result.support,
1260 }),
1261 error: None,
1262 }),
1263 Err(err) => comparisons.push(MethodComparisonEntry {
1264 method: ScientificMethod::Eht,
1265 status: MethodComparisonStatus::Error,
1266 available: true,
1267 confidence: meta.confidence,
1268 confidence_score: meta.confidence_score,
1269 warnings: meta.warnings,
1270 limitations: meta.limitations,
1271 payload: None,
1272 error: Some(err),
1273 }),
1274 }
1275 }
1276 }
1277
1278 MethodComparisonResult { plan, comparisons }
1279}
1280
1281pub fn compute_esp(
1283 elements: &[u8],
1284 positions: &[[f64; 3]],
1285 spacing: f64,
1286 padding: f64,
1287) -> Result<esp::EspGrid, String> {
1288 let pop = compute_population(elements, positions)?;
1289 let (grid, _report) = gpu::esp_grid_gpu::compute_esp_grid_with_report(
1290 elements,
1291 positions,
1292 &pop.mulliken_charges,
1293 spacing,
1294 padding,
1295 );
1296 Ok(grid)
1297}
1298
1299pub fn compute_dos(
1301 elements: &[u8],
1302 positions: &[[f64; 3]],
1303 sigma: f64,
1304 e_min: f64,
1305 e_max: f64,
1306 n_points: usize,
1307) -> Result<dos::DosResult, String> {
1308 let eht_result = eht::solve_eht(elements, positions, None)?;
1309 let flat: Vec<f64> = positions.iter().flat_map(|p| p.iter().copied()).collect();
1310
1311 #[cfg(feature = "parallel")]
1312 {
1313 Ok(dos::compute_pdos_parallel(
1314 elements,
1315 &flat,
1316 &eht_result.energies,
1317 &eht_result.coefficients,
1318 eht_result.n_electrons,
1319 sigma,
1320 e_min,
1321 e_max,
1322 n_points,
1323 ))
1324 }
1325
1326 #[cfg(not(feature = "parallel"))]
1327 {
1328 Ok(dos::compute_pdos(
1329 elements,
1330 &flat,
1331 &eht_result.energies,
1332 &eht_result.coefficients,
1333 eht_result.n_electrons,
1334 sigma,
1335 e_min,
1336 e_max,
1337 n_points,
1338 ))
1339 }
1340}
1341
1342pub fn compute_rmsd(coords: &[f64], reference: &[f64]) -> f64 {
1344 alignment::compute_rmsd(coords, reference)
1345}
1346
1347pub fn compute_uff_energy(smiles: &str, coords: &[f64]) -> Result<f64, String> {
1349 let mol = graph::Molecule::from_smiles(smiles)?;
1350 let n = mol.graph.node_count();
1351 if coords.len() != n * 3 {
1352 return Err(format!("coords length {} != 3 * atoms {}", coords.len(), n));
1353 }
1354 let ff = forcefield::builder::build_uff_force_field(&mol);
1355 let mut gradient = vec![0.0f64; n * 3];
1356 let energy = ff.compute_system_energy_and_gradients(coords, &mut gradient);
1357 Ok(energy)
1358}
1359
1360pub fn compute_mmff94_energy(smiles: &str, coords: &[f64]) -> Result<f64, String> {
1367 let mol = graph::Molecule::from_smiles(smiles)?;
1368 let n = mol.graph.node_count();
1369 if coords.len() != n * 3 {
1370 return Err(format!("coords length {} != 3 * atoms {}", coords.len(), n));
1371 }
1372 let elements: Vec<u8> = (0..n)
1373 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
1374 .collect();
1375 let bonds: Vec<(usize, usize, u8)> = mol
1376 .graph
1377 .edge_indices()
1378 .map(|e| {
1379 let (a, b) = mol.graph.edge_endpoints(e).unwrap();
1380 let order = match mol.graph[e].order {
1381 graph::BondOrder::Single => 1u8,
1382 graph::BondOrder::Double => 2,
1383 graph::BondOrder::Triple => 3,
1384 graph::BondOrder::Aromatic => 2,
1385 graph::BondOrder::Unknown => 1,
1386 };
1387 (a.index(), b.index(), order)
1388 })
1389 .collect();
1390 let terms = forcefield::mmff94::Mmff94Builder::build(&elements, &bonds);
1391 let (energy, _grad) = forcefield::mmff94::Mmff94Builder::total_energy(&terms, coords);
1392 Ok(energy)
1393}
1394
1395pub fn compute_pm3(elements: &[u8], positions: &[[f64; 3]]) -> Result<pm3::Pm3Result, String> {
1402 pm3::solve_pm3(elements, positions)
1403}
1404
1405pub fn compute_pm3_gradient(
1410 elements: &[u8],
1411 positions: &[[f64; 3]],
1412) -> Result<pm3::Pm3GradientResult, String> {
1413 pm3::compute_pm3_gradient(elements, positions)
1414}
1415
1416pub fn compute_xtb(elements: &[u8], positions: &[[f64; 3]]) -> Result<xtb::XtbResult, String> {
1423 xtb::solve_xtb(elements, positions)
1424}
1425
1426pub fn compute_xtb_gradient(
1431 elements: &[u8],
1432 positions: &[[f64; 3]],
1433) -> Result<xtb::XtbGradientResult, String> {
1434 xtb::compute_xtb_gradient(elements, positions)
1435}
1436
1437pub fn compute_stda_uvvis(
1444 elements: &[u8],
1445 positions: &[[f64; 3]],
1446 sigma: f64,
1447 e_min: f64,
1448 e_max: f64,
1449 n_points: usize,
1450 broadening: reactivity::BroadeningType,
1451) -> Result<reactivity::StdaUvVisSpectrum, String> {
1452 reactivity::compute_stda_uvvis_spectrum(
1453 elements, positions, sigma, e_min, e_max, n_points, broadening,
1454 )
1455}
1456
1457pub fn compute_vibrational_analysis(
1464 elements: &[u8],
1465 positions: &[[f64; 3]],
1466 method: &str,
1467 step_size: Option<f64>,
1468) -> Result<ir::VibrationalAnalysis, String> {
1469 let hessian_method =
1470 match method {
1471 "eht" => ir::HessianMethod::Eht,
1472 "pm3" => ir::HessianMethod::Pm3,
1473 "xtb" => ir::HessianMethod::Xtb,
1474 "uff" => return Err(
1475 "UFF vibrational analysis requires SMILES; use compute_vibrational_analysis_uff"
1476 .to_string(),
1477 ),
1478 _ => return Err(format!("Unknown method '{}', use eht/pm3/xtb/uff", method)),
1479 };
1480 ir::compute_vibrational_analysis(elements, positions, hessian_method, step_size)
1481}
1482
1483pub fn compute_vibrational_analysis_uff(
1488 smiles: &str,
1489 elements: &[u8],
1490 positions: &[[f64; 3]],
1491 step_size: Option<f64>,
1492) -> Result<ir::VibrationalAnalysis, String> {
1493 ir::vibrations::compute_vibrational_analysis_uff(smiles, elements, positions, step_size)
1494}
1495
1496pub fn compute_ir_spectrum(
1502 analysis: &ir::VibrationalAnalysis,
1503 gamma: f64,
1504 wn_min: f64,
1505 wn_max: f64,
1506 n_points: usize,
1507) -> ir::IrSpectrum {
1508 ir::compute_ir_spectrum(analysis, gamma, wn_min, wn_max, n_points)
1509}
1510
1511pub fn predict_nmr_shifts(smiles: &str) -> Result<nmr::NmrShiftResult, String> {
1513 let mol = graph::Molecule::from_smiles(smiles)?;
1514 Ok(nmr::predict_chemical_shifts(&mol))
1515}
1516
1517pub fn predict_nmr_shifts_for_nucleus(
1519 smiles: &str,
1520 nucleus: &str,
1521) -> Result<Vec<nmr::ChemicalShift>, String> {
1522 let mol = graph::Molecule::from_smiles(smiles)?;
1523 let nucleus = parse_nmr_nucleus(nucleus)?;
1524 Ok(nmr::predict_chemical_shifts_for_nucleus(&mol, nucleus))
1525}
1526
1527pub fn predict_nmr_couplings(
1532 smiles: &str,
1533 positions: &[[f64; 3]],
1534) -> Result<Vec<nmr::JCoupling>, String> {
1535 let mol = graph::Molecule::from_smiles(smiles)?;
1536 Ok(nmr::predict_j_couplings(&mol, positions))
1537}
1538
1539fn parse_nmr_nucleus(nucleus: &str) -> Result<nmr::NmrNucleus, String> {
1540 nmr::NmrNucleus::parse(nucleus)
1541}
1542
1543struct GiaoNmrPreflight {
1544 nucleus: nmr::NmrNucleus,
1545 system: scf::types::MolecularSystem,
1546 basis: scf::basis::BasisSet,
1547 fallback_elements: Vec<u8>,
1548}
1549
1550const EXPLICIT_GIAO_STO3G_ELEMENTS: &[u8] =
1551 &[1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 14, 15, 16, 17, 35, 53];
1552
1553fn has_explicit_giao_sto3g_support(z: u8) -> bool {
1554 EXPLICIT_GIAO_STO3G_ELEMENTS.contains(&z)
1555}
1556
1557fn format_element_label(z: u8) -> String {
1558 nmr::NmrNucleus::default_for_element(z)
1559 .map(|nucleus| format!("{}({})", nucleus.element_symbol(), z))
1560 .unwrap_or_else(|| format!("Z{}", z))
1561}
1562
1563fn format_element_set(elements: &[u8]) -> String {
1564 let mut unique = elements.to_vec();
1565 unique.sort_unstable();
1566 unique.dedup();
1567 unique
1568 .into_iter()
1569 .map(format_element_label)
1570 .collect::<Vec<_>>()
1571 .join(", ")
1572}
1573
1574fn preflight_giao_nmr(
1575 elements: &[u8],
1576 positions: &[[f64; 3]],
1577 nucleus: &str,
1578 config: &GiaoNmrConfig,
1579) -> Result<GiaoNmrPreflight, String> {
1580 if elements.is_empty() {
1581 return Err("GIAO NMR requires at least one atom.".to_string());
1582 }
1583 if elements.len() != positions.len() {
1584 return Err(format!(
1585 "GIAO NMR requires one 3D coordinate per atom; got {} elements and {} positions.",
1586 elements.len(),
1587 positions.len()
1588 ));
1589 }
1590 if config.max_scf_iterations == 0 {
1591 return Err("GIAO NMR requires max_scf_iterations > 0.".to_string());
1592 }
1593 if config.multiplicity != 1 {
1594 return Err(format!(
1595 "The public GIAO NMR API currently supports singlet closed-shell systems only; got multiplicity {}.",
1596 config.multiplicity
1597 ));
1598 }
1599
1600 let parsed_nucleus = parse_nmr_nucleus(nucleus)?;
1601 if !elements
1602 .iter()
1603 .any(|&z| z == parsed_nucleus.atomic_number())
1604 {
1605 return Err(format!(
1606 "Requested nucleus {} is not present in the provided element list.",
1607 parsed_nucleus.canonical()
1608 ));
1609 }
1610
1611 let system = scf::types::MolecularSystem::from_angstrom(
1612 elements,
1613 positions,
1614 config.charge,
1615 config.multiplicity,
1616 );
1617 if !system.n_electrons().is_multiple_of(2) {
1618 return Err(format!(
1619 "The public GIAO NMR API currently supports even-electron closed-shell systems only; got {} electrons.",
1620 system.n_electrons()
1621 ));
1622 }
1623
1624 let fallback_elements: Vec<u8> = elements
1625 .iter()
1626 .copied()
1627 .filter(|&z| !has_explicit_giao_sto3g_support(z))
1628 .collect();
1629
1630 if !config.allow_basis_fallback && !fallback_elements.is_empty() {
1631 return Err(format!(
1632 "The public GIAO NMR SCF route only has explicit STO-3G basis data for {}. This system includes fallback-only elements: {}. Re-run with GiaoNmrConfig.allow_basis_fallback=true for a screening-only attempt.",
1633 format_element_set(EXPLICIT_GIAO_STO3G_ELEMENTS),
1634 format_element_set(&fallback_elements)
1635 ));
1636 }
1637
1638 let basis = scf::basis::BasisSet::sto3g(&system.atomic_numbers, &system.positions_bohr);
1639 let n_occupied = system.n_electrons() / 2;
1640 if n_occupied > basis.n_basis {
1641 return Err(format!(
1642 "GIAO NMR preflight failed for {}: STO-3G supplies {} basis functions for {} occupied orbitals. This usually means one or more elements only have the hydrogen-like fallback basis in the current SCF path.",
1643 parsed_nucleus.canonical(),
1644 basis.n_basis,
1645 n_occupied
1646 ));
1647 }
1648
1649 Ok(GiaoNmrPreflight {
1650 nucleus: parsed_nucleus,
1651 system,
1652 basis,
1653 fallback_elements,
1654 })
1655}
1656
1657pub fn compute_giao_nmr(
1659 elements: &[u8],
1660 positions: &[[f64; 3]],
1661 nucleus: &str,
1662) -> Result<GiaoNmrResult, String> {
1663 compute_giao_nmr_configured(elements, positions, nucleus, &GiaoNmrConfig::default())
1664}
1665
1666pub fn compute_giao_nmr_configured(
1668 elements: &[u8],
1669 positions: &[[f64; 3]],
1670 nucleus: &str,
1671 config: &GiaoNmrConfig,
1672) -> Result<GiaoNmrResult, String> {
1673 let request = preflight_giao_nmr(elements, positions, nucleus, config)?;
1674
1675 let scf_config = scf::scf_loop::ScfConfig {
1676 max_iterations: config.max_scf_iterations,
1677 use_parallel_eri: config.use_parallel_eri,
1678 parallel_threshold: if config.use_parallel_eri { 0 } else { 20 },
1679 ..scf::scf_loop::ScfConfig::default()
1680 };
1681
1682 let scf = scf::scf_loop::run_scf(&request.system, &scf_config);
1683 let shieldings = spectroscopy::compute_nmr_shieldings_for_nucleus(
1684 request.nucleus,
1685 &request.system.atomic_numbers,
1686 &request.system.positions_bohr,
1687 &spectroscopy::ScfInput::from(&scf),
1688 &request.basis.function_to_atom,
1689 );
1690 if shieldings.is_empty() {
1691 return Err(format!(
1692 "Requested nucleus {} is not present in the provided element list.",
1693 request.nucleus.canonical()
1694 ));
1695 }
1696
1697 let entries: Vec<GiaoNmrEntry> = shieldings
1698 .iter()
1699 .map(|entry| GiaoNmrEntry {
1700 atom_index: entry.atom_index,
1701 element: entry.element,
1702 tensor: entry.tensor,
1703 isotropic: entry.isotropic,
1704 anisotropy: entry.anisotropy,
1705 chemical_shift: entry.chemical_shift,
1706 })
1707 .collect();
1708
1709 let mut notes = vec![format!(
1710 "GIAO NMR for {} via the public SCF route (RHF/STO-3G).",
1711 request.nucleus.unicode_label()
1712 )];
1713 if scf.converged {
1714 notes.push(format!(
1715 "SCF converged in {} iterations with total energy {:.6} Hartree.",
1716 scf.scf_iterations, scf.total_energy
1717 ));
1718 } else {
1719 notes.push(format!(
1720 "SCF did not reach the requested threshold after {} iterations; treat the returned shieldings as screening-level only.",
1721 scf.scf_iterations
1722 ));
1723 }
1724 if request.fallback_elements.is_empty() {
1725 notes.push(
1726 "All elements in this system use explicit STO-3G basis data in the current SCF path."
1727 .to_string(),
1728 );
1729 } else {
1730 notes.push(format!(
1731 "Fallback basis enabled for {}. Heavy-element shieldings are qualitative in this mode.",
1732 format_element_set(&request.fallback_elements)
1733 ));
1734 }
1735 if request.nucleus.is_quadrupolar() {
1736 notes.push(
1737 "Quadrupolar nucleus selected: this API returns isotropic shieldings and chemical shifts only; relaxation-driven broadening still needs experimental lineshape modeling.".to_string(),
1738 );
1739 }
1740 if !request.nucleus.is_primary_target() {
1741 notes.push(
1742 "Non-1H/13C nuclei remain sensitive to the reference shielding and basis choice; benchmark the target family before quantitative use.".to_string(),
1743 );
1744 }
1745
1746 Ok(GiaoNmrResult {
1747 nucleus: request.nucleus.canonical().to_string(),
1748 target_atomic_number: request.nucleus.atomic_number(),
1749 method: "RHF/GIAO".to_string(),
1750 basis_set: "STO-3G".to_string(),
1751 charge: config.charge,
1752 multiplicity: config.multiplicity,
1753 scf_converged: scf.converged,
1754 scf_iterations: scf.scf_iterations,
1755 total_energy_hartree: scf.total_energy,
1756 n_basis: scf.n_basis,
1757 n_target_atoms: entries.len(),
1758 chemical_shifts: entries.iter().map(|entry| entry.chemical_shift).collect(),
1759 shieldings: entries,
1760 fallback_elements: request.fallback_elements,
1761 notes,
1762 })
1763}
1764
1765pub fn compute_nmr_spectrum(
1772 smiles: &str,
1773 nucleus: &str,
1774 gamma: f64,
1775 ppm_min: f64,
1776 ppm_max: f64,
1777 n_points: usize,
1778) -> Result<nmr::NmrSpectrum, String> {
1779 compute_nmr_spectrum_with_coords(smiles, &[], nucleus, gamma, ppm_min, ppm_max, n_points)
1780}
1781
1782pub fn compute_nmr_spectrum_with_coords(
1785 smiles: &str,
1786 positions: &[[f64; 3]],
1787 nucleus: &str,
1788 gamma: f64,
1789 ppm_min: f64,
1790 ppm_max: f64,
1791 n_points: usize,
1792) -> Result<nmr::NmrSpectrum, String> {
1793 let mol = graph::Molecule::from_smiles(smiles)?;
1794 let couplings = nmr::predict_j_couplings(&mol, positions);
1795 let nuc = parse_nmr_nucleus(nucleus)?;
1796 let shifts = nmr::predict_chemical_shifts_for_nucleus(&mol, nuc);
1797 Ok(nmr::spectrum::compute_nmr_spectrum_for_shifts(
1798 &shifts, &couplings, nuc, gamma, ppm_min, ppm_max, n_points,
1799 ))
1800}
1801
1802pub fn compute_hose_codes(smiles: &str, max_radius: usize) -> Result<Vec<nmr::HoseCode>, String> {
1804 let mol = graph::Molecule::from_smiles(smiles)?;
1805 Ok(nmr::hose::generate_hose_codes(&mol, max_radius))
1806}
1807
1808pub fn compute_orbital_mesh(
1820 elements: &[u8],
1821 positions: &[[f64; 3]],
1822 method: &str,
1823 mo_index: usize,
1824 spacing: f64,
1825 padding: f64,
1826 isovalue: f32,
1827) -> Result<mesh::OrbitalMeshResult, String> {
1828 let m = match method.to_lowercase().as_str() {
1829 "eht" | "huckel" => mesh::MeshMethod::Eht,
1830 "pm3" => mesh::MeshMethod::Pm3,
1831 "xtb" | "gfn0" | "gfn-xtb" => mesh::MeshMethod::Xtb,
1832 "hf3c" | "hf-3c" | "hf" => mesh::MeshMethod::Hf3c,
1833 _ => {
1834 return Err(format!(
1835 "Unknown method '{}'. Supported: eht, pm3, xtb, hf3c",
1836 method
1837 ))
1838 }
1839 };
1840 mesh::compute_orbital_mesh(elements, positions, m, mo_index, spacing, padding, isovalue)
1841}
1842
1843pub fn compute_ml_descriptors(
1850 elements: &[u8],
1851 bonds: &[(usize, usize, u8)],
1852 charges: &[f64],
1853 aromatic_atoms: &[bool],
1854) -> ml::MolecularDescriptors {
1855 ml::compute_descriptors(elements, bonds, charges, aromatic_atoms)
1856}
1857
1858pub fn predict_ml_properties(desc: &ml::MolecularDescriptors) -> ml::MlPropertyResult {
1863 ml::predict_properties(desc)
1864}
1865
1866pub fn predict_ensemble(
1871 elements: &[u8],
1872 bonds: &[(usize, usize, u8)],
1873 charges: &[f64],
1874 aromatic_atoms: &[bool],
1875) -> ml::EnsembleResult {
1876 let desc = ml::compute_descriptors(elements, bonds, charges, aromatic_atoms);
1877 ml::predict_ensemble(&desc, elements, bonds)
1878}
1879
1880pub fn compute_tpsa(elements: &[u8], bonds: &[(usize, usize, u8)]) -> f64 {
1882 ml::compute_tpsa(elements, bonds)
1883}
1884
1885pub fn compute_md_trajectory(
1887 smiles: &str,
1888 coords: &[f64],
1889 n_steps: usize,
1890 dt_fs: f64,
1891 seed: u64,
1892) -> Result<dynamics::MdTrajectory, String> {
1893 dynamics::simulate_velocity_verlet_uff(smiles, coords, n_steps, dt_fs, seed, None)
1894}
1895
1896pub fn compute_md_trajectory_nvt(
1898 smiles: &str,
1899 coords: &[f64],
1900 n_steps: usize,
1901 dt_fs: f64,
1902 seed: u64,
1903 target_temp_k: f64,
1904 thermostat_tau_fs: f64,
1905) -> Result<dynamics::MdTrajectory, String> {
1906 dynamics::simulate_velocity_verlet_uff(
1907 smiles,
1908 coords,
1909 n_steps,
1910 dt_fs,
1911 seed,
1912 Some((target_temp_k, thermostat_tau_fs)),
1913 )
1914}
1915
1916pub fn compute_simplified_neb_path(
1918 smiles: &str,
1919 start_coords: &[f64],
1920 end_coords: &[f64],
1921 n_images: usize,
1922 n_iter: usize,
1923 spring_k: f64,
1924 step_size: f64,
1925) -> Result<dynamics::NebPathResult, String> {
1926 dynamics::compute_simplified_neb_path(
1927 smiles,
1928 start_coords,
1929 end_coords,
1930 n_images,
1931 n_iter,
1932 spring_k,
1933 step_size,
1934 )
1935}
1936
1937fn coords_flat_to_matrix_f32(coords: &[f64], n_atoms: usize) -> nalgebra::DMatrix<f32> {
1938 let mut m = nalgebra::DMatrix::<f32>::zeros(n_atoms, 3);
1939 for i in 0..n_atoms {
1940 m[(i, 0)] = coords[3 * i] as f32;
1941 m[(i, 1)] = coords[3 * i + 1] as f32;
1942 m[(i, 2)] = coords[3 * i + 2] as f32;
1943 }
1944 m
1945}
1946
1947fn coords_matrix_f32_to_flat(m: &nalgebra::DMatrix<f32>) -> Vec<f64> {
1948 let mut out = Vec::with_capacity(m.nrows() * 3);
1949 for i in 0..m.nrows() {
1950 out.push(m[(i, 0)] as f64);
1951 out.push(m[(i, 1)] as f64);
1952 out.push(m[(i, 2)] as f64);
1953 }
1954 out
1955}
1956
1957pub fn search_conformers_with_uff(
1966 smiles: &str,
1967 n_samples: usize,
1968 seed: u64,
1969 rmsd_threshold: f64,
1970) -> Result<ConformerSearchResult, String> {
1971 if n_samples == 0 {
1972 return Err("n_samples must be > 0".to_string());
1973 }
1974 if rmsd_threshold <= 0.0 {
1975 return Err("rmsd_threshold must be > 0".to_string());
1976 }
1977
1978 let mol = graph::Molecule::from_smiles(smiles)?;
1979 let n_atoms = mol.graph.node_count();
1980 let bounds = distgeom::smooth_bounds_matrix(distgeom::calculate_bounds_matrix(&mol));
1981
1982 let mut generated = Vec::new();
1983 let mut notes = Vec::new();
1984 let mut rotatable_bonds = 0usize;
1985
1986 for i in 0..n_samples {
1987 let sample_seed = seed.wrapping_add(i as u64);
1988 let conf = embed(smiles, sample_seed);
1989
1990 if conf.error.is_some() || conf.coords.len() != n_atoms * 3 {
1991 continue;
1992 }
1993
1994 let mut coords = coords_flat_to_matrix_f32(&conf.coords, n_atoms);
1995 let rot_mc = forcefield::optimize_torsions_monte_carlo_bounds(
1996 &mut coords,
1997 &mol,
1998 &bounds,
1999 sample_seed ^ 0x9E37_79B9_7F4A_7C15,
2000 64,
2001 0.4,
2002 );
2003 let rot_greedy = forcefield::optimize_torsions_bounds(&mut coords, &mol, &bounds, 2);
2004 let rot = rot_mc.max(rot_greedy);
2005 rotatable_bonds = rot;
2006 let coords_flat = coords_matrix_f32_to_flat(&coords);
2007
2008 match compute_uff_energy(smiles, &coords_flat) {
2009 Ok(energy_kcal_mol) => generated.push(ConformerEnsembleMember {
2010 seed: sample_seed,
2011 cluster_id: None,
2012 coords: coords_flat,
2013 energy_kcal_mol,
2014 }),
2015 Err(_) => continue,
2016 }
2017 }
2018
2019 if generated.is_empty() {
2020 return Err("failed to generate any valid conformers".to_string());
2021 }
2022
2023 generated.sort_by(|a, b| {
2024 a.energy_kcal_mol
2025 .partial_cmp(&b.energy_kcal_mol)
2026 .unwrap_or(std::cmp::Ordering::Equal)
2027 });
2028
2029 let generated_count = generated.len();
2030
2031 let mut unique = Vec::new();
2032 let mut cluster_members: Vec<Vec<u64>> = Vec::new();
2033 for candidate in generated {
2034 let existing_cluster = unique.iter().position(|u: &ConformerEnsembleMember| {
2035 compute_rmsd(&candidate.coords, &u.coords) < rmsd_threshold
2036 });
2037
2038 if let Some(cluster_id) = existing_cluster {
2039 cluster_members[cluster_id].push(candidate.seed);
2040 } else {
2041 unique.push(candidate.clone());
2042 cluster_members.push(vec![candidate.seed]);
2043 }
2044 }
2045
2046 for (cluster_id, member) in unique.iter_mut().enumerate() {
2047 member.cluster_id = Some(cluster_id);
2048 }
2049
2050 let clusters: Vec<ConformerClusterSummary> = unique
2051 .iter()
2052 .enumerate()
2053 .map(|(cluster_id, representative)| ConformerClusterSummary {
2054 cluster_id,
2055 representative_seed: representative.seed,
2056 size: cluster_members[cluster_id].len(),
2057 member_seeds: cluster_members[cluster_id].clone(),
2058 })
2059 .collect();
2060
2061 notes.push(
2062 "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."
2063 .to_string(),
2064 );
2065
2066 Ok(ConformerSearchResult {
2067 generated: generated_count,
2068 unique: unique.len(),
2069 rotatable_bonds,
2070 conformers: unique,
2071 clusters,
2072 notes,
2073 })
2074}
2075
2076pub fn compute_uff_energy_with_aromatic_heuristics(
2078 smiles: &str,
2079 coords: &[f64],
2080) -> Result<reactivity::UffHeuristicEnergy, String> {
2081 let mol = graph::Molecule::from_smiles(smiles)?;
2082 let n = mol.graph.node_count();
2083 if coords.len() != n * 3 {
2084 return Err(format!("coords length {} != 3 * atoms {}", coords.len(), n));
2085 }
2086
2087 let ff = forcefield::builder::build_uff_force_field(&mol);
2088 let mut gradient = vec![0.0f64; n * 3];
2089 let raw = ff.compute_system_energy_and_gradients(coords, &mut gradient);
2090 Ok(reactivity::apply_aromatic_uff_correction(&mol, raw))
2091}
2092
2093pub fn compute_empirical_pka(smiles: &str) -> Result<reactivity::EmpiricalPkaResult, String> {
2095 let mol = graph::Molecule::from_smiles(smiles)?;
2096 let charges = compute_charges(smiles)?;
2097 Ok(reactivity::estimate_empirical_pka(&mol, &charges.charges))
2098}
2099
2100pub fn create_unit_cell(
2102 a: f64,
2103 b: f64,
2104 c: f64,
2105 alpha: f64,
2106 beta: f64,
2107 gamma: f64,
2108) -> materials::UnitCell {
2109 materials::UnitCell::from_parameters(&materials::CellParameters {
2110 a,
2111 b,
2112 c,
2113 alpha,
2114 beta,
2115 gamma,
2116 })
2117}
2118
2119pub fn assemble_framework(
2123 node: &materials::Sbu,
2124 linker: &materials::Sbu,
2125 topology: &materials::Topology,
2126 cell: &materials::UnitCell,
2127) -> materials::CrystalStructure {
2128 materials::assemble_framework(node, linker, topology, cell)
2129}
2130
2131pub fn compute_hf3c(
2139 elements: &[u8],
2140 positions: &[[f64; 3]],
2141 config: &hf::HfConfig,
2142) -> Result<hf::Hf3cResult, String> {
2143 hf::solve_hf3c(elements, positions, config)
2144}
2145
2146pub fn compute_ani(elements: &[u8], positions: &[[f64; 3]]) -> Result<ani::AniResult, String> {
2154 ani::api::compute_ani_test(elements, positions)
2155}
2156
2157pub fn compute_ani_with_models(
2164 elements: &[u8],
2165 positions: &[[f64; 3]],
2166 config: &ani::AniConfig,
2167 models: &std::collections::HashMap<u8, ani::nn::FeedForwardNet>,
2168) -> Result<ani::AniResult, String> {
2169 ani::compute_ani(elements, positions, config, models)
2170}
2171
2172pub fn compute_esp_grid(
2176 elements: &[u8],
2177 positions: &[[f64; 3]],
2178 spacing: f64,
2179 padding: f64,
2180) -> Result<esp::EspGrid, String> {
2181 compute_esp(elements, positions, spacing, padding)
2182}
2183
2184pub fn to_canonical_smiles(smiles: &str) -> Result<String, String> {
2191 let mol = graph::Molecule::from_smiles(smiles)?;
2192 Ok(canonical_smiles::to_canonical_smiles(&mol))
2193}
2194
2195pub fn analyze_stereo(smiles: &str, coords: &[f64]) -> Result<stereo::StereoAnalysis, String> {
2199 let mol = graph::Molecule::from_smiles(smiles)?;
2200 let positions: Vec<[f64; 3]> = coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
2201 Ok(stereo::analyze_stereo(&mol, &positions))
2202}
2203
2204pub fn compute_nonpolar_solvation(
2208 elements: &[u8],
2209 positions: &[[f64; 3]],
2210 probe_radius: Option<f64>,
2211) -> solvation::NonPolarSolvation {
2212 solvation::compute_nonpolar_solvation(elements, positions, probe_radius)
2213}
2214
2215pub fn compute_gb_solvation(
2217 elements: &[u8],
2218 positions: &[[f64; 3]],
2219 charges: &[f64],
2220 solvent_dielectric: Option<f64>,
2221 solute_dielectric: Option<f64>,
2222 probe_radius: Option<f64>,
2223) -> solvation::GbSolvation {
2224 solvation::compute_gb_solvation(
2225 elements,
2226 positions,
2227 charges,
2228 solvent_dielectric,
2229 solute_dielectric,
2230 probe_radius,
2231 )
2232}
2233
2234pub fn butina_cluster(conformers: &[Vec<f64>], rmsd_cutoff: f64) -> clustering::ClusterResult {
2238 clustering::butina_cluster(conformers, rmsd_cutoff)
2239}
2240
2241pub fn compute_rmsd_matrix(conformers: &[Vec<f64>]) -> Vec<Vec<f64>> {
2243 clustering::compute_rmsd_matrix(conformers)
2244}
2245
2246pub fn compute_sssr(smiles: &str) -> Result<rings::sssr::SssrResult, String> {
2250 let mol = graph::Molecule::from_smiles(smiles)?;
2251 Ok(rings::sssr::compute_sssr(&mol))
2252}
2253
2254pub fn compute_ecfp(
2256 smiles: &str,
2257 radius: usize,
2258 n_bits: usize,
2259) -> Result<rings::ecfp::ECFPFingerprint, String> {
2260 let mol = graph::Molecule::from_smiles(smiles)?;
2261 Ok(rings::ecfp::compute_ecfp(&mol, radius, n_bits))
2262}
2263
2264pub fn compute_ecfp_batch(
2266 smiles_list: &[&str],
2267 radius: usize,
2268 n_bits: usize,
2269) -> Result<Vec<rings::ecfp::ECFPFingerprint>, String> {
2270 let mols: Result<Vec<_>, _> = smiles_list
2271 .iter()
2272 .map(|s| graph::Molecule::from_smiles(s))
2273 .collect();
2274 let mols = mols?;
2275 Ok(rings::ecfp::compute_ecfp_batch(&mols, radius, n_bits))
2276}
2277
2278pub fn compute_tanimoto(
2280 fp1: &rings::ecfp::ECFPFingerprint,
2281 fp2: &rings::ecfp::ECFPFingerprint,
2282) -> f64 {
2283 rings::ecfp::compute_tanimoto(fp1, fp2)
2284}
2285
2286pub fn compute_charges_configured(
2290 smiles: &str,
2291 config: &charges::gasteiger::GasteigerConfig,
2292) -> Result<charges::gasteiger::ChargeResult, String> {
2293 let mol = graph::Molecule::from_smiles(smiles)?;
2294 let n = mol.graph.node_count();
2295 let elements: Vec<u8> = (0..n)
2296 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
2297 .collect();
2298 let bonds: Vec<(usize, usize)> = mol
2299 .graph
2300 .edge_indices()
2301 .map(|e| {
2302 let (a, b) = mol.graph.edge_endpoints(e).unwrap();
2303 (a.index(), b.index())
2304 })
2305 .collect();
2306 let formal_charges: Vec<i8> = (0..n)
2307 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].formal_charge)
2308 .collect();
2309 charges::gasteiger::gasteiger_marsili_charges_configured(
2310 &elements,
2311 &bonds,
2312 &formal_charges,
2313 config,
2314 )
2315}
2316
2317pub fn compute_ensemble_j_couplings(
2321 smiles: &str,
2322 conformer_coords: &[Vec<f64>],
2323 energies_kcal: &[f64],
2324 temperature_k: f64,
2325) -> Result<Vec<nmr::JCoupling>, String> {
2326 let mol = graph::Molecule::from_smiles(smiles)?;
2327 let positions: Vec<Vec<[f64; 3]>> = conformer_coords
2328 .iter()
2329 .map(|c| c.chunks(3).map(|p| [p[0], p[1], p[2]]).collect())
2330 .collect();
2331 Ok(nmr::coupling::ensemble_averaged_j_couplings(
2332 &mol,
2333 &positions,
2334 energies_kcal,
2335 temperature_k,
2336 ))
2337}
2338
2339pub fn compute_ir_spectrum_broadened(
2343 analysis: &ir::VibrationalAnalysis,
2344 gamma: f64,
2345 wn_min: f64,
2346 wn_max: f64,
2347 n_points: usize,
2348 broadening: &str,
2349) -> ir::IrSpectrum {
2350 let bt = match broadening {
2351 "gaussian" | "Gaussian" => ir::BroadeningType::Gaussian,
2352 _ => ir::BroadeningType::Lorentzian,
2353 };
2354 ir::vibrations::compute_ir_spectrum_with_broadening(
2355 analysis, gamma, wn_min, wn_max, n_points, bt,
2356 )
2357}
2358
2359#[cfg(test)]
2360mod tests {
2361 use super::*;
2362
2363 #[test]
2364 fn test_cisplatin_style_support_metadata() {
2365 let smiles = "[Pt](Cl)(Cl)([NH3])[NH3]";
2366 let mol = parse(smiles).expect("Cisplatin-style example should parse");
2367 let elements: Vec<u8> = (0..mol.graph.node_count())
2368 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
2369 .collect();
2370
2371 let caps = get_system_capabilities(&elements);
2372 assert!(caps.eht.available);
2373 assert!(matches!(
2374 caps.eht.confidence,
2375 eht::SupportLevel::Experimental
2376 ));
2377 }
2378
2379 #[test]
2380 fn test_pt_system_routes_to_uff_when_experimental_disabled() {
2381 let smiles = "[Pt](Cl)(Cl)([NH3])[NH3]";
2382 let mol = parse(smiles).expect("Pt example should parse");
2383 let n = mol.graph.node_count();
2384 let elements: Vec<u8> = (0..n)
2385 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
2386 .collect();
2387 let positions = vec![[0.0, 0.0, 0.0]; n];
2388
2389 let result = compute_eht_or_uff_fallback(smiles, &elements, &positions, false)
2390 .expect("Fallback workflow should return a result");
2391
2392 assert!(matches!(
2393 result,
2394 ElectronicWorkflowResult::UffFallback { .. }
2395 ));
2396 }
2397
2398 #[test]
2399 fn test_method_plan_prefers_supported_workflows_for_organic_systems() {
2400 let plan = get_system_method_plan(&[6, 1, 1, 1, 1]);
2401
2402 assert_eq!(plan.geometry.recommended, Some(ScientificMethod::Embed));
2403 assert_eq!(
2404 plan.force_field_energy.recommended,
2405 Some(ScientificMethod::Uff)
2406 );
2407 assert_eq!(plan.orbitals.recommended, Some(ScientificMethod::Eht));
2408 assert_eq!(plan.population.recommended, Some(ScientificMethod::Eht));
2409 assert_eq!(plan.orbital_grid.recommended, Some(ScientificMethod::Eht));
2410 assert!(plan.orbitals.methods[0].confidence_score > 0.9);
2411 }
2412
2413 #[test]
2414 fn test_method_plan_marks_metal_orbital_workflow_experimental() {
2415 let plan = get_system_method_plan(&[78, 17, 17, 7, 7]);
2416
2417 assert_eq!(plan.orbitals.recommended, Some(ScientificMethod::Eht));
2418 assert_eq!(
2419 plan.force_field_energy.recommended,
2420 Some(ScientificMethod::Uff)
2421 );
2422 assert!(matches!(
2423 plan.orbitals.methods[0].confidence,
2424 eht::SupportLevel::Experimental
2425 ));
2426 assert!(plan.orbitals.methods[0].confidence_score < 0.9);
2427 assert!(!plan.orbitals.methods[0].warnings.is_empty());
2428 }
2429
2430 #[test]
2431 fn test_method_plan_reports_unavailable_workflows_for_unsupported_elements() {
2432 let plan = get_system_method_plan(&[92]);
2433
2434 assert_eq!(plan.force_field_energy.recommended, None);
2435 assert_eq!(plan.orbitals.recommended, None);
2436 assert_eq!(plan.population.recommended, None);
2437 assert_eq!(plan.orbital_grid.recommended, None);
2438 assert!(!plan.orbitals.methods[0].limitations.is_empty());
2439 }
2440
2441 #[test]
2442 fn test_compare_methods_supported_system_returns_success_rows() {
2443 let result = compare_methods("CC", &[6, 6], &[[0.0, 0.0, 0.0], [1.54, 0.0, 0.0]], true);
2444 assert_eq!(result.comparisons.len(), 2);
2445 assert!(result
2446 .comparisons
2447 .iter()
2448 .any(|entry| matches!(entry.method, ScientificMethod::Uff) && entry.available));
2449 assert!(result.comparisons.iter().any(|entry| matches!(
2450 entry.method,
2451 ScientificMethod::Eht
2452 ) && matches!(
2453 entry.status,
2454 MethodComparisonStatus::Success
2455 )));
2456 }
2457
2458 #[test]
2459 fn test_compare_methods_blocks_experimental_eht_when_disabled() {
2460 let result = compare_methods("[O]", &[78], &[[0.0, 0.0, 0.0]], false);
2461 let eht_row = result
2462 .comparisons
2463 .iter()
2464 .find(|entry| matches!(entry.method, ScientificMethod::Eht))
2465 .expect("EHT row must exist");
2466 assert!(matches!(
2467 eht_row.status,
2468 MethodComparisonStatus::Unavailable
2469 ));
2470 }
2471
2472 #[test]
2473 fn test_compare_methods_reports_unavailable_for_unsupported_elements() {
2474 let result = compare_methods("[O]", &[92], &[[0.0, 0.0, 0.0]], true);
2475 assert!(result
2476 .comparisons
2477 .iter()
2478 .all(|entry| matches!(entry.status, MethodComparisonStatus::Unavailable)));
2479 }
2480
2481 #[test]
2482 fn test_compute_fukui_descriptors_returns_atomwise_output() {
2483 let elements = [8u8, 1, 1];
2484 let positions = [[0.0, 0.0, 0.0], [0.757, 0.586, 0.0], [-0.757, 0.586, 0.0]];
2485 let result = compute_fukui_descriptors(&elements, &positions).unwrap();
2486 assert_eq!(result.num_atoms, 3);
2487 assert_eq!(result.condensed.len(), 3);
2488 }
2489
2490 #[test]
2491 fn test_compute_uv_vis_spectrum_returns_requested_grid_size() {
2492 let elements = [6u8, 6, 1, 1, 1, 1];
2493 let positions = [
2494 [0.0, 0.0, 0.0],
2495 [1.34, 0.0, 0.0],
2496 [-0.6, 0.92, 0.0],
2497 [-0.6, -0.92, 0.0],
2498 [1.94, 0.92, 0.0],
2499 [1.94, -0.92, 0.0],
2500 ];
2501 let spectrum = compute_uv_vis_spectrum(&elements, &positions, 0.2, 0.5, 8.0, 256).unwrap();
2502 assert_eq!(spectrum.energies_ev.len(), 256);
2503 assert_eq!(spectrum.intensities.len(), 256);
2504 }
2505
2506 #[test]
2507 fn test_analyze_graph_features_reports_benzene_aromaticity() {
2508 let analysis = analyze_graph_features("c1ccccc1").unwrap();
2509 assert_eq!(analysis.aromaticity.aromatic_atoms.len(), 12);
2510 assert_eq!(
2511 analysis
2512 .aromaticity
2513 .aromatic_atoms
2514 .iter()
2515 .filter(|v| **v)
2516 .count(),
2517 6
2518 );
2519 assert_eq!(analysis.aromaticity.aromatic_bonds.len(), 6);
2520 }
2521
2522 #[test]
2523 fn test_compute_empirical_pka_finds_acidic_site_for_acetic_acid() {
2524 let result = compute_empirical_pka("CC(=O)O").unwrap();
2525 assert!(!result.acidic_sites.is_empty());
2526 }
2527
2528 #[test]
2529 fn test_compute_uff_energy_with_aromatic_heuristics_applies_correction() {
2530 let conf = embed("c1ccccc1", 42);
2531 assert!(conf.error.is_none());
2532
2533 let result = compute_uff_energy_with_aromatic_heuristics("c1ccccc1", &conf.coords).unwrap();
2534 assert!(result.aromatic_bond_count >= 6);
2535 assert!(result.corrected_energy_kcal_mol <= result.raw_energy_kcal_mol);
2536 }
2537
2538 #[test]
2539 fn test_search_conformers_with_uff_returns_ranked_unique_ensemble() {
2540 let result = search_conformers_with_uff("CCCC", 10, 42, 0.2).unwrap();
2541 assert!(result.generated >= 1);
2542 assert!(result.unique >= 1);
2543 assert_eq!(result.unique, result.conformers.len());
2544 assert_eq!(result.unique, result.clusters.len());
2545 assert!(result.rotatable_bonds >= 1);
2546
2547 let mut total_members = 0usize;
2548 for (i, cluster) in result.clusters.iter().enumerate() {
2549 assert_eq!(cluster.cluster_id, i);
2550 assert!(cluster.size >= 1);
2551 total_members += cluster.size;
2552 assert_eq!(result.conformers[i].cluster_id, Some(i));
2553 assert_eq!(result.conformers[i].seed, cluster.representative_seed);
2554 }
2555 assert_eq!(total_members, result.generated);
2556
2557 for i in 1..result.conformers.len() {
2558 assert!(
2559 result.conformers[i - 1].energy_kcal_mol <= result.conformers[i].energy_kcal_mol
2560 );
2561 }
2562 }
2563
2564 #[test]
2565 fn test_search_conformers_with_uff_large_rmsd_threshold_collapses_duplicates() {
2566 let result = search_conformers_with_uff("CCCC", 8, 123, 10.0).unwrap();
2567 assert_eq!(result.unique, 1);
2568 assert_eq!(result.clusters.len(), 1);
2569 assert_eq!(result.clusters[0].size, result.generated);
2570 }
2571
2572 #[test]
2573 fn test_compute_md_trajectory_velocity_verlet_runs() {
2574 let conf = embed("CC", 42);
2575 assert!(conf.error.is_none());
2576
2577 let trj = compute_md_trajectory("CC", &conf.coords, 10, 0.25, 7).unwrap();
2578 assert_eq!(trj.frames.len(), 11);
2579 assert!(trj
2580 .frames
2581 .iter()
2582 .all(|f| f.coords.iter().all(|v| v.is_finite())));
2583 }
2584
2585 #[test]
2586 fn test_compute_md_trajectory_nvt_runs() {
2587 let conf = embed("CCO", 42);
2588 assert!(conf.error.is_none());
2589
2590 let trj =
2591 compute_md_trajectory_nvt("CCO", &conf.coords, 12, 0.25, 17, 300.0, 10.0).unwrap();
2592 assert_eq!(trj.frames.len(), 13);
2593 assert!(trj.frames.iter().all(|f| f.temperature_k.is_finite()));
2594 }
2595
2596 #[test]
2597 fn test_compute_simplified_neb_path_runs() {
2598 let c1 = embed("CC", 42);
2599 let c2 = embed("CC", 43);
2600 assert!(c1.error.is_none());
2601 assert!(c2.error.is_none());
2602
2603 let path =
2604 compute_simplified_neb_path("CC", &c1.coords, &c2.coords, 6, 20, 0.01, 1e-5).unwrap();
2605 assert_eq!(path.images.len(), 6);
2606 assert!(path
2607 .images
2608 .iter()
2609 .all(|img| img.potential_energy_kcal_mol.is_finite()));
2610 }
2611
2612 #[test]
2613 fn test_compute_hf3c_water() {
2614 let conf = embed("O", 42);
2615 assert!(conf.error.is_none());
2616 let pos: Vec<[f64; 3]> = conf.coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
2617
2618 let result = compute_hf3c(&conf.elements, &pos, &hf::HfConfig::default());
2619 assert!(result.is_ok(), "HF-3c should succeed for water");
2620 let r = result.unwrap();
2621 assert!(r.energy.is_finite());
2622 assert!(!r.orbital_energies.is_empty());
2623 }
2624
2625 #[test]
2626 fn test_compute_ani_water() {
2627 let conf = embed("O", 42);
2628 assert!(conf.error.is_none());
2629 let pos: Vec<[f64; 3]> = conf.coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
2630
2631 let result = compute_ani(&conf.elements, &pos);
2632 assert!(result.is_ok(), "ANI should succeed for water");
2633 let r = result.unwrap();
2634 assert!(r.energy.is_finite());
2635 assert_eq!(r.forces.len(), 3); assert_eq!(r.species.len(), 3);
2637 }
2638
2639 #[test]
2640 fn test_compute_esp_grid_water() {
2641 let conf = embed("O", 42);
2642 assert!(conf.error.is_none());
2643 let pos: Vec<[f64; 3]> = conf.coords.chunks(3).map(|c| [c[0], c[1], c[2]]).collect();
2644
2645 let result = compute_esp_grid(&conf.elements, &pos, 0.5, 3.0);
2646 assert!(result.is_ok(), "ESP grid should succeed for water");
2647 let g = result.unwrap();
2648 assert!(!g.values.is_empty());
2649 assert!(g.spacing > 0.0);
2650 assert!(g.dims[0] > 0 && g.dims[1] > 0 && g.dims[2] > 0);
2651 }
2652}