1#![allow(clippy::too_many_arguments)]
3#![allow(clippy::needless_range_loop)]
4
5pub mod alignment;
6pub mod charges;
7pub mod conformer;
8pub mod dipole;
9pub mod distgeom;
10pub mod dos;
11pub mod dynamics;
12pub mod eht;
13pub mod esp;
14pub mod etkdg;
15pub mod forcefield;
16pub mod graph;
17pub mod ir;
18pub mod materials;
19pub mod ml;
20pub mod nmr;
21pub mod optimization;
22pub mod pm3;
23pub mod population;
24pub mod reactivity;
25pub mod smarts;
26pub mod smiles;
27pub mod surface;
28pub mod topology;
29pub mod transport;
30pub mod xtb;
31
32use serde::{Deserialize, Serialize};
33
34#[derive(Debug, Clone, Serialize, Deserialize)]
38pub struct ConformerResult {
39 pub smiles: String,
41 pub num_atoms: usize,
43 pub coords: Vec<f64>,
46 pub elements: Vec<u8>,
48 pub bonds: Vec<(usize, usize, String)>,
50 pub error: Option<String>,
52 pub time_ms: f64,
54}
55
56#[derive(Debug, Clone, Serialize, Deserialize)]
58pub struct ConformerConfig {
59 pub seed: u64,
61 pub num_threads: usize,
63}
64
65#[derive(Debug, Clone, Serialize, Deserialize)]
67pub struct ConformerEnsembleMember {
68 pub seed: u64,
70 pub cluster_id: Option<usize>,
72 pub coords: Vec<f64>,
74 pub energy_kcal_mol: f64,
76}
77
78#[derive(Debug, Clone, Serialize, Deserialize)]
80pub struct ConformerClusterSummary {
81 pub cluster_id: usize,
83 pub representative_seed: u64,
85 pub size: usize,
87 pub member_seeds: Vec<u64>,
89}
90
91#[derive(Debug, Clone, Serialize, Deserialize)]
93pub struct ConformerSearchResult {
94 pub generated: usize,
96 pub unique: usize,
98 pub rotatable_bonds: usize,
100 pub conformers: Vec<ConformerEnsembleMember>,
102 pub clusters: Vec<ConformerClusterSummary>,
104 pub notes: Vec<String>,
106}
107
108#[derive(Debug, Clone, Serialize, Deserialize)]
110pub struct MethodCapability {
111 pub available: bool,
113 pub confidence: eht::SupportLevel,
115 pub unsupported_elements: Vec<u8>,
117 pub warnings: Vec<String>,
119}
120
121#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
123#[serde(rename_all = "snake_case")]
124pub enum ScientificMethod {
125 Embed,
126 Uff,
127 Eht,
128 Pm3,
129 Xtb,
130 Mmff94,
131}
132
133#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
135#[serde(rename_all = "snake_case")]
136pub enum PropertyRequest {
137 Geometry,
138 ForceFieldEnergy,
139 Orbitals,
140 Population,
141 OrbitalGrid,
142}
143
144#[derive(Debug, Clone, Serialize, Deserialize)]
146pub struct MethodMetadata {
147 pub method: ScientificMethod,
148 pub available: bool,
149 pub confidence: eht::SupportLevel,
150 pub confidence_score: f64,
151 pub limitations: Vec<String>,
152 pub warnings: Vec<String>,
153}
154
155#[derive(Debug, Clone, Serialize, Deserialize)]
157pub struct PropertyMethodPlan {
158 pub property: PropertyRequest,
159 pub recommended: Option<ScientificMethod>,
160 pub fallback: Option<ScientificMethod>,
161 pub rationale: Vec<String>,
162 pub methods: Vec<MethodMetadata>,
163}
164
165#[derive(Debug, Clone, Serialize, Deserialize)]
167pub struct SystemMethodPlan {
168 pub capabilities: SystemCapabilities,
169 pub geometry: PropertyMethodPlan,
170 pub force_field_energy: PropertyMethodPlan,
171 pub orbitals: PropertyMethodPlan,
172 pub population: PropertyMethodPlan,
173 pub orbital_grid: PropertyMethodPlan,
174}
175
176#[derive(Debug, Clone, Serialize, Deserialize)]
178pub struct SystemCapabilities {
179 pub embed: MethodCapability,
180 pub uff: MethodCapability,
181 pub eht: MethodCapability,
182 pub population: MethodCapability,
183 pub orbital_grid: MethodCapability,
184}
185
186#[derive(Debug, Clone, Serialize, Deserialize)]
188#[serde(tag = "mode", rename_all = "snake_case")]
189pub enum ElectronicWorkflowResult {
190 Eht {
191 result: eht::EhtResult,
192 },
193 UffFallback {
194 energy_kcal_mol: f64,
195 reason: String,
196 support: eht::EhtSupport,
197 },
198}
199
200#[derive(Debug, Clone, Copy, PartialEq, Eq, Serialize, Deserialize)]
202#[serde(rename_all = "snake_case")]
203pub enum MethodComparisonStatus {
204 Success,
205 Unavailable,
206 Error,
207}
208
209#[derive(Debug, Clone, Serialize, Deserialize)]
211#[serde(tag = "kind", rename_all = "snake_case")]
212pub enum MethodComparisonPayload {
213 Eht {
214 homo_energy: f64,
215 lumo_energy: f64,
216 gap: f64,
217 support: eht::EhtSupport,
218 },
219 Uff {
220 energy_kcal_mol: f64,
221 },
222}
223
224#[derive(Debug, Clone, Serialize, Deserialize)]
226pub struct MethodComparisonEntry {
227 pub method: ScientificMethod,
228 pub status: MethodComparisonStatus,
229 pub available: bool,
230 pub confidence: eht::SupportLevel,
231 pub confidence_score: f64,
232 pub warnings: Vec<String>,
233 pub limitations: Vec<String>,
234 pub payload: Option<MethodComparisonPayload>,
235 pub error: Option<String>,
236}
237
238#[derive(Debug, Clone, Serialize, Deserialize)]
240pub struct MethodComparisonResult {
241 pub plan: SystemMethodPlan,
242 pub comparisons: Vec<MethodComparisonEntry>,
243}
244
245#[derive(Debug, Clone, Serialize, Deserialize)]
247pub struct AromaticityAnalysis {
248 pub aromatic_atoms: Vec<bool>,
249 pub aromatic_bonds: Vec<(usize, usize)>,
250}
251
252#[derive(Debug, Clone, Serialize, Deserialize)]
254pub struct StereocenterAnalysis {
255 pub tagged_tetrahedral_centers: Vec<usize>,
256 pub inferred_tetrahedral_centers: Vec<usize>,
257}
258
259#[derive(Debug, Clone, Serialize, Deserialize)]
261pub struct GraphFeatureAnalysis {
262 pub aromaticity: AromaticityAnalysis,
263 pub stereocenters: StereocenterAnalysis,
264}
265
266impl Default for ConformerConfig {
267 fn default() -> Self {
268 Self {
269 seed: 42,
270 num_threads: 0,
271 }
272 }
273}
274
275pub fn version() -> String {
279 format!("sci-form {}", env!("CARGO_PKG_VERSION"))
280}
281
282pub fn get_eht_support(elements: &[u8]) -> eht::EhtSupport {
284 eht::analyze_eht_support(elements)
285}
286
287fn is_uff_element_supported(z: u8) -> bool {
288 matches!(
289 z,
290 1 | 5
291 | 6
292 | 7
293 | 8
294 | 9
295 | 14
296 | 15
297 | 16
298 | 17
299 | 22
300 | 23
301 | 24
302 | 25
303 | 26
304 | 27
305 | 28
306 | 29
307 | 30
308 | 32
309 | 33
310 | 34
311 | 35
312 | 42
313 | 46
314 | 47
315 | 50
316 | 51
317 | 52
318 | 53
319 | 78
320 | 79
321 )
322}
323
324fn unique_sorted_unsupported(elements: &[u8], pred: impl Fn(u8) -> bool) -> Vec<u8> {
325 let mut out: Vec<u8> = elements.iter().copied().filter(|&z| !pred(z)).collect();
326 out.sort_unstable();
327 out.dedup();
328 out
329}
330
331pub fn get_system_capabilities(elements: &[u8]) -> SystemCapabilities {
333 let eht_support = get_eht_support(elements);
334 let uff_unsupported = unique_sorted_unsupported(elements, is_uff_element_supported);
335
336 let embed = MethodCapability {
337 available: !elements.is_empty(),
338 confidence: eht::SupportLevel::Experimental,
339 unsupported_elements: Vec::new(),
340 warnings: vec![
341 "Embed capability is inferred from element presence only; final success still depends on full molecular graph and geometry constraints.".to_string(),
342 ],
343 };
344
345 let uff = if uff_unsupported.is_empty() {
346 MethodCapability {
347 available: true,
348 confidence: eht::SupportLevel::Supported,
349 unsupported_elements: Vec::new(),
350 warnings: Vec::new(),
351 }
352 } else {
353 MethodCapability {
354 available: false,
355 confidence: eht::SupportLevel::Unsupported,
356 unsupported_elements: uff_unsupported.clone(),
357 warnings: vec![format!(
358 "UFF atom typing is unavailable for elements {:?}.",
359 uff_unsupported
360 )],
361 }
362 };
363
364 let eht = MethodCapability {
365 available: eht_support.unsupported_elements.is_empty(),
366 confidence: eht_support.level,
367 unsupported_elements: eht_support.unsupported_elements.clone(),
368 warnings: eht_support.warnings.clone(),
369 };
370
371 let population = MethodCapability {
372 available: eht.available,
373 confidence: eht.confidence,
374 unsupported_elements: eht.unsupported_elements.clone(),
375 warnings: eht.warnings.clone(),
376 };
377
378 let orbital_grid = MethodCapability {
379 available: eht.available,
380 confidence: eht.confidence,
381 unsupported_elements: eht.unsupported_elements.clone(),
382 warnings: eht.warnings.clone(),
383 };
384
385 SystemCapabilities {
386 embed,
387 uff,
388 eht,
389 population,
390 orbital_grid,
391 }
392}
393
394fn confidence_score_for_method(method: ScientificMethod, capability: &MethodCapability) -> f64 {
395 if !capability.available {
396 return 0.0;
397 }
398
399 match method {
400 ScientificMethod::Embed => 0.8,
401 ScientificMethod::Uff | ScientificMethod::Mmff94 => match capability.confidence {
402 eht::SupportLevel::Supported => 0.95,
403 eht::SupportLevel::Experimental => 0.75,
404 eht::SupportLevel::Unsupported => 0.0,
405 },
406 ScientificMethod::Eht | ScientificMethod::Pm3 | ScientificMethod::Xtb => {
407 match capability.confidence {
408 eht::SupportLevel::Supported => 0.95,
409 eht::SupportLevel::Experimental => 0.6,
410 eht::SupportLevel::Unsupported => 0.0,
411 }
412 }
413 }
414}
415
416fn build_method_metadata(
417 method: ScientificMethod,
418 capability: &MethodCapability,
419 extra_limitations: &[&str],
420) -> MethodMetadata {
421 let mut limitations: Vec<String> = extra_limitations.iter().map(|s| s.to_string()).collect();
422
423 if !capability.unsupported_elements.is_empty() {
424 limitations.push(format!(
425 "Unsupported elements for this method: {:?}.",
426 capability.unsupported_elements
427 ));
428 }
429
430 if matches!(method, ScientificMethod::Eht)
431 && matches!(capability.confidence, eht::SupportLevel::Experimental)
432 {
433 limitations.push(
434 "Transition-metal EHT parameters remain provisional and should be treated as experimental."
435 .to_string(),
436 );
437 }
438
439 MethodMetadata {
440 method,
441 available: capability.available,
442 confidence: capability.confidence,
443 confidence_score: confidence_score_for_method(method, capability),
444 limitations,
445 warnings: capability.warnings.clone(),
446 }
447}
448
449fn build_property_plan(
450 property: PropertyRequest,
451 recommended: Option<ScientificMethod>,
452 fallback: Option<ScientificMethod>,
453 rationale: Vec<String>,
454 methods: Vec<MethodMetadata>,
455) -> PropertyMethodPlan {
456 PropertyMethodPlan {
457 property,
458 recommended,
459 fallback,
460 rationale,
461 methods,
462 }
463}
464
465pub fn get_system_method_plan(elements: &[u8]) -> SystemMethodPlan {
467 let capabilities = get_system_capabilities(elements);
468
469 let geometry_method = build_method_metadata(
470 ScientificMethod::Embed,
471 &capabilities.embed,
472 &["Geometry generation still depends on graph topology, stereochemistry, and embedding constraints."],
473 );
474 let geometry = build_property_plan(
475 PropertyRequest::Geometry,
476 geometry_method.available.then_some(ScientificMethod::Embed),
477 None,
478 vec!["Embedding is the top-level geometry generation path in sci-form.".to_string()],
479 vec![geometry_method],
480 );
481
482 let uff_method = build_method_metadata(
483 ScientificMethod::Uff,
484 &capabilities.uff,
485 &["This recommendation applies to force-field energy evaluation, not molecular orbital analysis."],
486 );
487 let force_field_energy = build_property_plan(
488 PropertyRequest::ForceFieldEnergy,
489 uff_method.available.then_some(ScientificMethod::Uff),
490 None,
491 vec![
492 "UFF is the top-level force-field energy path when atom typing is available."
493 .to_string(),
494 ],
495 vec![uff_method],
496 );
497
498 let eht_method = build_method_metadata(
499 ScientificMethod::Eht,
500 &capabilities.eht,
501 &["EHT is the only current top-level orbital method in sci-form."],
502 );
503 let orbitals = build_property_plan(
504 PropertyRequest::Orbitals,
505 eht_method.available.then_some(ScientificMethod::Eht),
506 None,
507 vec!["Orbital energies and MO coefficients are produced by the EHT workflow.".to_string()],
508 vec![eht_method.clone()],
509 );
510
511 let population_method = build_method_metadata(
512 ScientificMethod::Eht,
513 &capabilities.population,
514 &["Population analysis is derived from the EHT density and overlap matrices."],
515 );
516 let population = build_property_plan(
517 PropertyRequest::Population,
518 population_method.available.then_some(ScientificMethod::Eht),
519 None,
520 vec!["Population analysis currently requires a successful EHT calculation.".to_string()],
521 vec![population_method],
522 );
523
524 let orbital_grid_method = build_method_metadata(
525 ScientificMethod::Eht,
526 &capabilities.orbital_grid,
527 &["Orbital-grid rendering currently depends on EHT molecular-orbital coefficients."],
528 );
529 let orbital_grid = build_property_plan(
530 PropertyRequest::OrbitalGrid,
531 orbital_grid_method
532 .available
533 .then_some(ScientificMethod::Eht),
534 None,
535 vec![
536 "Orbital-grid generation currently requires a successful EHT calculation.".to_string(),
537 ],
538 vec![orbital_grid_method],
539 );
540
541 SystemMethodPlan {
542 capabilities,
543 geometry,
544 force_field_energy,
545 orbitals,
546 population,
547 orbital_grid,
548 }
549}
550
551pub fn embed(smiles: &str, seed: u64) -> ConformerResult {
553 #[cfg(not(target_arch = "wasm32"))]
554 let start = std::time::Instant::now();
555
556 let mol = match graph::Molecule::from_smiles(smiles) {
557 Ok(m) => m,
558 Err(e) => {
559 return ConformerResult {
560 smiles: smiles.to_string(),
561 num_atoms: 0,
562 coords: vec![],
563 elements: vec![],
564 bonds: vec![],
565 error: Some(e),
566 #[cfg(not(target_arch = "wasm32"))]
567 time_ms: start.elapsed().as_secs_f64() * 1000.0,
568 #[cfg(target_arch = "wasm32")]
569 time_ms: 0.0,
570 };
571 }
572 };
573
574 let n = mol.graph.node_count();
575 let elements: Vec<u8> = (0..n)
576 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
577 .collect();
578 let bonds: Vec<(usize, usize, String)> = mol
579 .graph
580 .edge_indices()
581 .map(|e| {
582 let (a, b) = mol.graph.edge_endpoints(e).unwrap();
583 let order = match mol.graph[e].order {
584 graph::BondOrder::Single => "SINGLE",
585 graph::BondOrder::Double => "DOUBLE",
586 graph::BondOrder::Triple => "TRIPLE",
587 graph::BondOrder::Aromatic => "AROMATIC",
588 graph::BondOrder::Unknown => "UNKNOWN",
589 };
590 (a.index(), b.index(), order.to_string())
591 })
592 .collect();
593
594 match conformer::generate_3d_conformer(&mol, seed) {
595 Ok(coords) => {
596 let mut flat = Vec::with_capacity(n * 3);
597 for i in 0..n {
598 flat.push(coords[(i, 0)] as f64);
599 flat.push(coords[(i, 1)] as f64);
600 flat.push(coords[(i, 2)] as f64);
601 }
602 ConformerResult {
603 smiles: smiles.to_string(),
604 num_atoms: n,
605 coords: flat,
606 elements,
607 bonds,
608 error: None,
609 #[cfg(not(target_arch = "wasm32"))]
610 time_ms: start.elapsed().as_secs_f64() * 1000.0,
611 #[cfg(target_arch = "wasm32")]
612 time_ms: 0.0,
613 }
614 }
615 Err(e) => ConformerResult {
616 smiles: smiles.to_string(),
617 num_atoms: n,
618 coords: vec![],
619 elements,
620 bonds,
621 error: Some(e),
622 #[cfg(not(target_arch = "wasm32"))]
623 time_ms: start.elapsed().as_secs_f64() * 1000.0,
624 #[cfg(target_arch = "wasm32")]
625 time_ms: 0.0,
626 },
627 }
628}
629
630#[cfg(feature = "parallel")]
635pub fn embed_batch(smiles_list: &[&str], config: &ConformerConfig) -> Vec<ConformerResult> {
636 use rayon::prelude::*;
637
638 if config.num_threads > 0 {
639 let pool = rayon::ThreadPoolBuilder::new()
640 .num_threads(config.num_threads)
641 .build()
642 .unwrap();
643 pool.install(|| {
644 smiles_list
645 .par_iter()
646 .map(|smi| embed(smi, config.seed))
647 .collect()
648 })
649 } else {
650 smiles_list
651 .par_iter()
652 .map(|smi| embed(smi, config.seed))
653 .collect()
654 }
655}
656
657#[cfg(not(feature = "parallel"))]
659pub fn embed_batch(smiles_list: &[&str], config: &ConformerConfig) -> Vec<ConformerResult> {
660 smiles_list
661 .iter()
662 .map(|smi| embed(smi, config.seed))
663 .collect()
664}
665
666pub fn parse(smiles: &str) -> Result<graph::Molecule, String> {
668 graph::Molecule::from_smiles(smiles)
669}
670
671pub fn compute_charges(smiles: &str) -> Result<charges::gasteiger::ChargeResult, String> {
676 let mol = graph::Molecule::from_smiles(smiles)?;
677 let n = mol.graph.node_count();
678 let elements: Vec<u8> = (0..n)
679 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
680 .collect();
681 let formal_charges: Vec<i8> = (0..n)
682 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].formal_charge)
683 .collect();
684 let bonds: Vec<(usize, usize)> = mol
685 .graph
686 .edge_indices()
687 .map(|e| {
688 let (a, b) = mol.graph.edge_endpoints(e).unwrap();
689 (a.index(), b.index())
690 })
691 .collect();
692 charges::gasteiger::gasteiger_marsili_charges(&elements, &bonds, &formal_charges, 6)
693}
694
695pub fn compute_sasa(
699 elements: &[u8],
700 coords_flat: &[f64],
701 probe_radius: Option<f64>,
702) -> Result<surface::sasa::SasaResult, String> {
703 if coords_flat.len() != elements.len() * 3 {
704 return Err(format!(
705 "coords length {} != 3 * elements {}",
706 coords_flat.len(),
707 elements.len()
708 ));
709 }
710 let positions: Vec<[f64; 3]> = coords_flat
711 .chunks_exact(3)
712 .map(|c| [c[0], c[1], c[2]])
713 .collect();
714
715 #[cfg(feature = "parallel")]
716 {
717 Ok(surface::sasa::compute_sasa_parallel(
718 elements,
719 &positions,
720 probe_radius,
721 None,
722 ))
723 }
724
725 #[cfg(not(feature = "parallel"))]
726 {
727 Ok(surface::sasa::compute_sasa(
728 elements,
729 &positions,
730 probe_radius,
731 None,
732 ))
733 }
734}
735
736pub fn compute_population(
738 elements: &[u8],
739 positions: &[[f64; 3]],
740) -> Result<population::PopulationResult, String> {
741 let eht_result = eht::solve_eht(elements, positions, None)?;
742 Ok(population::compute_population(
743 elements,
744 positions,
745 &eht_result.coefficients,
746 eht_result.n_electrons,
747 ))
748}
749
750pub fn compute_dipole(
752 elements: &[u8],
753 positions: &[[f64; 3]],
754) -> Result<dipole::DipoleResult, String> {
755 let eht_result = eht::solve_eht(elements, positions, None)?;
756 Ok(dipole::compute_dipole_from_eht(
757 elements,
758 positions,
759 &eht_result.coefficients,
760 eht_result.n_electrons,
761 ))
762}
763
764pub fn compute_frontier_descriptors(
766 elements: &[u8],
767 positions: &[[f64; 3]],
768) -> Result<reactivity::FrontierDescriptors, String> {
769 let eht_result = eht::solve_eht(elements, positions, None)?;
770 Ok(reactivity::compute_frontier_descriptors(
771 elements,
772 positions,
773 &eht_result,
774 ))
775}
776
777pub fn compute_fukui_descriptors(
779 elements: &[u8],
780 positions: &[[f64; 3]],
781) -> Result<reactivity::FukuiDescriptors, String> {
782 let eht_result = eht::solve_eht(elements, positions, None)?;
783 Ok(reactivity::compute_fukui_descriptors(
784 elements,
785 positions,
786 &eht_result,
787 ))
788}
789
790pub fn compute_reactivity_ranking(
792 elements: &[u8],
793 positions: &[[f64; 3]],
794) -> Result<reactivity::ReactivityRanking, String> {
795 let eht_result = eht::solve_eht(elements, positions, None)?;
796 let fukui = reactivity::compute_fukui_descriptors(elements, positions, &eht_result);
797 let pop = population::compute_population(
798 elements,
799 positions,
800 &eht_result.coefficients,
801 eht_result.n_electrons,
802 );
803 Ok(reactivity::rank_reactivity_sites(
804 &fukui,
805 &pop.mulliken_charges,
806 ))
807}
808
809pub fn compute_uv_vis_spectrum(
811 elements: &[u8],
812 positions: &[[f64; 3]],
813 sigma: f64,
814 e_min: f64,
815 e_max: f64,
816 n_points: usize,
817) -> Result<reactivity::UvVisSpectrum, String> {
818 let eht_result = eht::solve_eht(elements, positions, None)?;
819 Ok(reactivity::compute_uv_vis_like_spectrum(
820 &eht_result,
821 sigma,
822 e_min,
823 e_max,
824 n_points,
825 ))
826}
827
828pub fn compute_bond_orders(
830 elements: &[u8],
831 positions: &[[f64; 3]],
832) -> Result<population::BondOrderResult, String> {
833 let eht_result = eht::solve_eht(elements, positions, None)?;
834 Ok(population::compute_bond_orders(
835 elements,
836 positions,
837 &eht_result.coefficients,
838 eht_result.n_electrons,
839 ))
840}
841
842pub fn compute_topology(
844 elements: &[u8],
845 positions: &[[f64; 3]],
846) -> topology::TopologyAnalysisResult {
847 topology::analyze_topology(elements, positions)
848}
849
850pub fn analyze_graph_features(smiles: &str) -> Result<GraphFeatureAnalysis, String> {
852 use petgraph::visit::EdgeRef;
853
854 let mol = parse(smiles)?;
855 let n_atoms = mol.graph.node_count();
856 let mut aromatic_atoms = vec![false; n_atoms];
857 let mut aromatic_bonds = Vec::new();
858
859 for edge in mol.graph.edge_references() {
860 if matches!(edge.weight().order, graph::BondOrder::Aromatic) {
861 let i = edge.source().index();
862 let j = edge.target().index();
863 aromatic_atoms[i] = true;
864 aromatic_atoms[j] = true;
865 aromatic_bonds.push((i, j));
866 }
867 }
868
869 let mut tagged_tetrahedral_centers = Vec::new();
870 let mut inferred_tetrahedral_centers = Vec::new();
871 for i in 0..n_atoms {
872 let idx = petgraph::graph::NodeIndex::new(i);
873 let atom = &mol.graph[idx];
874 if matches!(
875 atom.chiral_tag,
876 graph::ChiralType::TetrahedralCW | graph::ChiralType::TetrahedralCCW
877 ) {
878 tagged_tetrahedral_centers.push(i);
879 }
880
881 let neighbors: Vec<_> = mol.graph.neighbors(idx).collect();
882 if neighbors.len() == 4 && matches!(atom.hybridization, graph::Hybridization::SP3) {
883 let mut signature: Vec<u8> = neighbors.iter().map(|n| mol.graph[*n].element).collect();
884 signature.sort_unstable();
885 signature.dedup();
886 if signature.len() >= 3 {
887 inferred_tetrahedral_centers.push(i);
888 }
889 }
890 }
891
892 Ok(GraphFeatureAnalysis {
893 aromaticity: AromaticityAnalysis {
894 aromatic_atoms,
895 aromatic_bonds,
896 },
897 stereocenters: StereocenterAnalysis {
898 tagged_tetrahedral_centers,
899 inferred_tetrahedral_centers,
900 },
901 })
902}
903
904pub fn compute_eht_or_uff_fallback(
909 smiles: &str,
910 elements: &[u8],
911 positions: &[[f64; 3]],
912 allow_experimental_eht: bool,
913) -> Result<ElectronicWorkflowResult, String> {
914 let support = get_eht_support(elements);
915 let should_fallback = match support.level {
916 eht::SupportLevel::Unsupported => true,
917 eht::SupportLevel::Experimental => !allow_experimental_eht,
918 eht::SupportLevel::Supported => false,
919 };
920
921 if should_fallback {
922 let coords_flat: Vec<f64> = positions.iter().flat_map(|p| p.iter().copied()).collect();
923 let energy = compute_uff_energy(smiles, &coords_flat).map_err(|e| {
924 format!(
925 "EHT is not appropriate for this system (support: {:?}) and UFF fallback failed: {}",
926 support.level, e
927 )
928 })?;
929 return Ok(ElectronicWorkflowResult::UffFallback {
930 energy_kcal_mol: energy,
931 reason: if matches!(support.level, eht::SupportLevel::Unsupported) {
932 "EHT unsupported for one or more elements; routed to UFF-only workflow.".to_string()
933 } else {
934 "EHT confidence is experimental and experimental mode is disabled; routed to UFF-only workflow."
935 .to_string()
936 },
937 support,
938 });
939 }
940
941 let result = eht::solve_eht(elements, positions, None)?;
942 Ok(ElectronicWorkflowResult::Eht { result })
943}
944
945pub fn compare_methods(
950 smiles: &str,
951 elements: &[u8],
952 positions: &[[f64; 3]],
953 allow_experimental_eht: bool,
954) -> MethodComparisonResult {
955 let plan = get_system_method_plan(elements);
956 let mut comparisons = Vec::new();
957
958 let coords_flat: Vec<f64> = positions.iter().flat_map(|p| p.iter().copied()).collect();
959
960 {
961 let meta = build_method_metadata(
962 ScientificMethod::Uff,
963 &plan.capabilities.uff,
964 &["Comparison uses UFF force-field energy as the UFF observable."],
965 );
966 if !meta.available {
967 comparisons.push(MethodComparisonEntry {
968 method: ScientificMethod::Uff,
969 status: MethodComparisonStatus::Unavailable,
970 available: false,
971 confidence: meta.confidence,
972 confidence_score: meta.confidence_score,
973 warnings: meta.warnings,
974 limitations: meta.limitations,
975 payload: None,
976 error: Some("UFF is unavailable for this element set.".to_string()),
977 });
978 } else {
979 match compute_uff_energy(smiles, &coords_flat) {
980 Ok(energy) => comparisons.push(MethodComparisonEntry {
981 method: ScientificMethod::Uff,
982 status: MethodComparisonStatus::Success,
983 available: true,
984 confidence: meta.confidence,
985 confidence_score: meta.confidence_score,
986 warnings: meta.warnings,
987 limitations: meta.limitations,
988 payload: Some(MethodComparisonPayload::Uff {
989 energy_kcal_mol: energy,
990 }),
991 error: None,
992 }),
993 Err(err) => comparisons.push(MethodComparisonEntry {
994 method: ScientificMethod::Uff,
995 status: MethodComparisonStatus::Error,
996 available: true,
997 confidence: meta.confidence,
998 confidence_score: meta.confidence_score,
999 warnings: meta.warnings,
1000 limitations: meta.limitations,
1001 payload: None,
1002 error: Some(err),
1003 }),
1004 }
1005 }
1006 }
1007
1008 {
1009 let meta = build_method_metadata(
1010 ScientificMethod::Eht,
1011 &plan.capabilities.eht,
1012 &["Comparison uses frontier orbital energies and gap as the EHT observable."],
1013 );
1014
1015 if !meta.available {
1016 comparisons.push(MethodComparisonEntry {
1017 method: ScientificMethod::Eht,
1018 status: MethodComparisonStatus::Unavailable,
1019 available: false,
1020 confidence: meta.confidence,
1021 confidence_score: meta.confidence_score,
1022 warnings: meta.warnings,
1023 limitations: meta.limitations,
1024 payload: None,
1025 error: Some("EHT is unavailable for this element set.".to_string()),
1026 });
1027 } else if matches!(meta.confidence, eht::SupportLevel::Experimental)
1028 && !allow_experimental_eht
1029 {
1030 comparisons.push(MethodComparisonEntry {
1031 method: ScientificMethod::Eht,
1032 status: MethodComparisonStatus::Unavailable,
1033 available: true,
1034 confidence: meta.confidence,
1035 confidence_score: meta.confidence_score,
1036 warnings: meta.warnings,
1037 limitations: meta.limitations,
1038 payload: None,
1039 error: Some(
1040 "EHT confidence is experimental and allow_experimental_eht=false.".to_string(),
1041 ),
1042 });
1043 } else {
1044 match eht::solve_eht(elements, positions, None) {
1045 Ok(result) => comparisons.push(MethodComparisonEntry {
1046 method: ScientificMethod::Eht,
1047 status: MethodComparisonStatus::Success,
1048 available: true,
1049 confidence: meta.confidence,
1050 confidence_score: meta.confidence_score,
1051 warnings: meta.warnings,
1052 limitations: meta.limitations,
1053 payload: Some(MethodComparisonPayload::Eht {
1054 homo_energy: result.homo_energy,
1055 lumo_energy: result.lumo_energy,
1056 gap: result.gap,
1057 support: result.support,
1058 }),
1059 error: None,
1060 }),
1061 Err(err) => comparisons.push(MethodComparisonEntry {
1062 method: ScientificMethod::Eht,
1063 status: MethodComparisonStatus::Error,
1064 available: true,
1065 confidence: meta.confidence,
1066 confidence_score: meta.confidence_score,
1067 warnings: meta.warnings,
1068 limitations: meta.limitations,
1069 payload: None,
1070 error: Some(err),
1071 }),
1072 }
1073 }
1074 }
1075
1076 MethodComparisonResult { plan, comparisons }
1077}
1078
1079pub fn compute_esp(
1081 elements: &[u8],
1082 positions: &[[f64; 3]],
1083 spacing: f64,
1084 padding: f64,
1085) -> Result<esp::EspGrid, String> {
1086 let pop = compute_population(elements, positions)?;
1087 #[cfg(feature = "parallel")]
1088 {
1089 Ok(esp::compute_esp_grid_parallel(
1090 elements,
1091 positions,
1092 &pop.mulliken_charges,
1093 spacing,
1094 padding,
1095 ))
1096 }
1097
1098 #[cfg(not(feature = "parallel"))]
1099 {
1100 Ok(esp::compute_esp_grid(
1101 elements,
1102 positions,
1103 &pop.mulliken_charges,
1104 spacing,
1105 padding,
1106 ))
1107 }
1108}
1109
1110pub fn compute_dos(
1112 elements: &[u8],
1113 positions: &[[f64; 3]],
1114 sigma: f64,
1115 e_min: f64,
1116 e_max: f64,
1117 n_points: usize,
1118) -> Result<dos::DosResult, String> {
1119 let eht_result = eht::solve_eht(elements, positions, None)?;
1120 let flat: Vec<f64> = positions.iter().flat_map(|p| p.iter().copied()).collect();
1121
1122 #[cfg(feature = "parallel")]
1123 {
1124 Ok(dos::compute_pdos_parallel(
1125 elements,
1126 &flat,
1127 &eht_result.energies,
1128 &eht_result.coefficients,
1129 eht_result.n_electrons,
1130 sigma,
1131 e_min,
1132 e_max,
1133 n_points,
1134 ))
1135 }
1136
1137 #[cfg(not(feature = "parallel"))]
1138 {
1139 Ok(dos::compute_pdos(
1140 elements,
1141 &flat,
1142 &eht_result.energies,
1143 &eht_result.coefficients,
1144 eht_result.n_electrons,
1145 sigma,
1146 e_min,
1147 e_max,
1148 n_points,
1149 ))
1150 }
1151}
1152
1153pub fn compute_rmsd(coords: &[f64], reference: &[f64]) -> f64 {
1155 alignment::compute_rmsd(coords, reference)
1156}
1157
1158pub fn compute_uff_energy(smiles: &str, coords: &[f64]) -> Result<f64, String> {
1160 let mol = graph::Molecule::from_smiles(smiles)?;
1161 let n = mol.graph.node_count();
1162 if coords.len() != n * 3 {
1163 return Err(format!("coords length {} != 3 * atoms {}", coords.len(), n));
1164 }
1165 let ff = forcefield::builder::build_uff_force_field(&mol);
1166 let mut gradient = vec![0.0f64; n * 3];
1167 let energy = ff.compute_system_energy_and_gradients(coords, &mut gradient);
1168 Ok(energy)
1169}
1170
1171pub fn compute_mmff94_energy(smiles: &str, coords: &[f64]) -> Result<f64, String> {
1178 let mol = graph::Molecule::from_smiles(smiles)?;
1179 let n = mol.graph.node_count();
1180 if coords.len() != n * 3 {
1181 return Err(format!("coords length {} != 3 * atoms {}", coords.len(), n));
1182 }
1183 let elements: Vec<u8> = (0..n)
1184 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
1185 .collect();
1186 let bonds: Vec<(usize, usize, u8)> = mol
1187 .graph
1188 .edge_indices()
1189 .map(|e| {
1190 let (a, b) = mol.graph.edge_endpoints(e).unwrap();
1191 let order = match mol.graph[e].order {
1192 graph::BondOrder::Single => 1u8,
1193 graph::BondOrder::Double => 2,
1194 graph::BondOrder::Triple => 3,
1195 graph::BondOrder::Aromatic => 2,
1196 graph::BondOrder::Unknown => 1,
1197 };
1198 (a.index(), b.index(), order)
1199 })
1200 .collect();
1201 let terms = forcefield::mmff94::Mmff94Builder::build(&elements, &bonds);
1202 let (energy, _grad) = forcefield::mmff94::Mmff94Builder::total_energy(&terms, coords);
1203 Ok(energy)
1204}
1205
1206pub fn compute_pm3(elements: &[u8], positions: &[[f64; 3]]) -> Result<pm3::Pm3Result, String> {
1213 pm3::solve_pm3(elements, positions)
1214}
1215
1216pub fn compute_xtb(elements: &[u8], positions: &[[f64; 3]]) -> Result<xtb::XtbResult, String> {
1223 xtb::solve_xtb(elements, positions)
1224}
1225
1226pub fn compute_stda_uvvis(
1233 elements: &[u8],
1234 positions: &[[f64; 3]],
1235 sigma: f64,
1236 e_min: f64,
1237 e_max: f64,
1238 n_points: usize,
1239 broadening: reactivity::BroadeningType,
1240) -> Result<reactivity::StdaUvVisSpectrum, String> {
1241 reactivity::compute_stda_uvvis_spectrum(
1242 elements, positions, sigma, e_min, e_max, n_points, broadening,
1243 )
1244}
1245
1246pub fn compute_vibrational_analysis(
1253 elements: &[u8],
1254 positions: &[[f64; 3]],
1255 method: &str,
1256 step_size: Option<f64>,
1257) -> Result<ir::VibrationalAnalysis, String> {
1258 let hessian_method = match method {
1259 "eht" => ir::HessianMethod::Eht,
1260 "pm3" => ir::HessianMethod::Pm3,
1261 "xtb" => ir::HessianMethod::Xtb,
1262 _ => return Err(format!("Unknown method '{}', use eht/pm3/xtb", method)),
1263 };
1264 ir::compute_vibrational_analysis(elements, positions, hessian_method, step_size)
1265}
1266
1267pub fn compute_ir_spectrum(
1273 analysis: &ir::VibrationalAnalysis,
1274 gamma: f64,
1275 wn_min: f64,
1276 wn_max: f64,
1277 n_points: usize,
1278) -> ir::IrSpectrum {
1279 ir::compute_ir_spectrum(analysis, gamma, wn_min, wn_max, n_points)
1280}
1281
1282pub fn predict_nmr_shifts(smiles: &str) -> Result<nmr::NmrShiftResult, String> {
1284 let mol = graph::Molecule::from_smiles(smiles)?;
1285 Ok(nmr::predict_chemical_shifts(&mol))
1286}
1287
1288pub fn predict_nmr_couplings(
1293 smiles: &str,
1294 positions: &[[f64; 3]],
1295) -> Result<Vec<nmr::JCoupling>, String> {
1296 let mol = graph::Molecule::from_smiles(smiles)?;
1297 Ok(nmr::predict_j_couplings(&mol, positions))
1298}
1299
1300pub fn compute_nmr_spectrum(
1307 smiles: &str,
1308 nucleus: &str,
1309 gamma: f64,
1310 ppm_min: f64,
1311 ppm_max: f64,
1312 n_points: usize,
1313) -> Result<nmr::NmrSpectrum, String> {
1314 let mol = graph::Molecule::from_smiles(smiles)?;
1315 let shifts = nmr::predict_chemical_shifts(&mol);
1316 let couplings = nmr::predict_j_couplings(&mol, &[]);
1317 let nuc = match nucleus {
1318 "1H" | "H1" | "h1" | "1h" | "proton" => nmr::NmrNucleus::H1,
1319 "13C" | "C13" | "c13" | "13c" | "carbon" => nmr::NmrNucleus::C13,
1320 _ => return Err(format!("Unknown nucleus '{}', use '1H' or '13C'", nucleus)),
1321 };
1322 Ok(nmr::compute_nmr_spectrum(
1323 &shifts, &couplings, nuc, gamma, ppm_min, ppm_max, n_points,
1324 ))
1325}
1326
1327pub fn compute_hose_codes(smiles: &str, max_radius: usize) -> Result<Vec<nmr::HoseCode>, String> {
1329 let mol = graph::Molecule::from_smiles(smiles)?;
1330 Ok(nmr::hose::generate_hose_codes(&mol, max_radius))
1331}
1332
1333pub fn compute_ml_descriptors(
1340 elements: &[u8],
1341 bonds: &[(usize, usize, u8)],
1342 charges: &[f64],
1343 aromatic_atoms: &[bool],
1344) -> ml::MolecularDescriptors {
1345 ml::compute_descriptors(elements, bonds, charges, aromatic_atoms)
1346}
1347
1348pub fn predict_ml_properties(desc: &ml::MolecularDescriptors) -> ml::MlPropertyResult {
1353 ml::predict_properties(desc)
1354}
1355
1356pub fn compute_md_trajectory(
1358 smiles: &str,
1359 coords: &[f64],
1360 n_steps: usize,
1361 dt_fs: f64,
1362 seed: u64,
1363) -> Result<dynamics::MdTrajectory, String> {
1364 dynamics::simulate_velocity_verlet_uff(smiles, coords, n_steps, dt_fs, seed, None)
1365}
1366
1367pub fn compute_md_trajectory_nvt(
1369 smiles: &str,
1370 coords: &[f64],
1371 n_steps: usize,
1372 dt_fs: f64,
1373 seed: u64,
1374 target_temp_k: f64,
1375 thermostat_tau_fs: f64,
1376) -> Result<dynamics::MdTrajectory, String> {
1377 dynamics::simulate_velocity_verlet_uff(
1378 smiles,
1379 coords,
1380 n_steps,
1381 dt_fs,
1382 seed,
1383 Some((target_temp_k, thermostat_tau_fs)),
1384 )
1385}
1386
1387pub fn compute_simplified_neb_path(
1389 smiles: &str,
1390 start_coords: &[f64],
1391 end_coords: &[f64],
1392 n_images: usize,
1393 n_iter: usize,
1394 spring_k: f64,
1395 step_size: f64,
1396) -> Result<dynamics::NebPathResult, String> {
1397 dynamics::compute_simplified_neb_path(
1398 smiles,
1399 start_coords,
1400 end_coords,
1401 n_images,
1402 n_iter,
1403 spring_k,
1404 step_size,
1405 )
1406}
1407
1408fn coords_flat_to_matrix_f32(coords: &[f64], n_atoms: usize) -> nalgebra::DMatrix<f32> {
1409 let mut m = nalgebra::DMatrix::<f32>::zeros(n_atoms, 3);
1410 for i in 0..n_atoms {
1411 m[(i, 0)] = coords[3 * i] as f32;
1412 m[(i, 1)] = coords[3 * i + 1] as f32;
1413 m[(i, 2)] = coords[3 * i + 2] as f32;
1414 }
1415 m
1416}
1417
1418fn coords_matrix_f32_to_flat(m: &nalgebra::DMatrix<f32>) -> Vec<f64> {
1419 let mut out = Vec::with_capacity(m.nrows() * 3);
1420 for i in 0..m.nrows() {
1421 out.push(m[(i, 0)] as f64);
1422 out.push(m[(i, 1)] as f64);
1423 out.push(m[(i, 2)] as f64);
1424 }
1425 out
1426}
1427
1428pub fn search_conformers_with_uff(
1437 smiles: &str,
1438 n_samples: usize,
1439 seed: u64,
1440 rmsd_threshold: f64,
1441) -> Result<ConformerSearchResult, String> {
1442 if n_samples == 0 {
1443 return Err("n_samples must be > 0".to_string());
1444 }
1445 if rmsd_threshold <= 0.0 {
1446 return Err("rmsd_threshold must be > 0".to_string());
1447 }
1448
1449 let mol = graph::Molecule::from_smiles(smiles)?;
1450 let n_atoms = mol.graph.node_count();
1451 let bounds = distgeom::smooth_bounds_matrix(distgeom::calculate_bounds_matrix(&mol));
1452
1453 let mut generated = Vec::new();
1454 let mut notes = Vec::new();
1455 let mut rotatable_bonds = 0usize;
1456
1457 for i in 0..n_samples {
1458 let sample_seed = seed.wrapping_add(i as u64);
1459 let conf = embed(smiles, sample_seed);
1460
1461 if conf.error.is_some() || conf.coords.len() != n_atoms * 3 {
1462 continue;
1463 }
1464
1465 let mut coords = coords_flat_to_matrix_f32(&conf.coords, n_atoms);
1466 let rot_mc = forcefield::optimize_torsions_monte_carlo_bounds(
1467 &mut coords,
1468 &mol,
1469 &bounds,
1470 sample_seed ^ 0x9E37_79B9_7F4A_7C15,
1471 64,
1472 0.4,
1473 );
1474 let rot_greedy = forcefield::optimize_torsions_bounds(&mut coords, &mol, &bounds, 2);
1475 let rot = rot_mc.max(rot_greedy);
1476 rotatable_bonds = rot;
1477 let coords_flat = coords_matrix_f32_to_flat(&coords);
1478
1479 match compute_uff_energy(smiles, &coords_flat) {
1480 Ok(energy_kcal_mol) => generated.push(ConformerEnsembleMember {
1481 seed: sample_seed,
1482 cluster_id: None,
1483 coords: coords_flat,
1484 energy_kcal_mol,
1485 }),
1486 Err(_) => continue,
1487 }
1488 }
1489
1490 if generated.is_empty() {
1491 return Err("failed to generate any valid conformers".to_string());
1492 }
1493
1494 generated.sort_by(|a, b| {
1495 a.energy_kcal_mol
1496 .partial_cmp(&b.energy_kcal_mol)
1497 .unwrap_or(std::cmp::Ordering::Equal)
1498 });
1499
1500 let generated_count = generated.len();
1501
1502 let mut unique = Vec::new();
1503 let mut cluster_members: Vec<Vec<u64>> = Vec::new();
1504 for candidate in generated {
1505 let existing_cluster = unique.iter().position(|u: &ConformerEnsembleMember| {
1506 compute_rmsd(&candidate.coords, &u.coords) < rmsd_threshold
1507 });
1508
1509 if let Some(cluster_id) = existing_cluster {
1510 cluster_members[cluster_id].push(candidate.seed);
1511 } else {
1512 unique.push(candidate.clone());
1513 cluster_members.push(vec![candidate.seed]);
1514 }
1515 }
1516
1517 for (cluster_id, member) in unique.iter_mut().enumerate() {
1518 member.cluster_id = Some(cluster_id);
1519 }
1520
1521 let clusters: Vec<ConformerClusterSummary> = unique
1522 .iter()
1523 .enumerate()
1524 .map(|(cluster_id, representative)| ConformerClusterSummary {
1525 cluster_id,
1526 representative_seed: representative.seed,
1527 size: cluster_members[cluster_id].len(),
1528 member_seeds: cluster_members[cluster_id].clone(),
1529 })
1530 .collect();
1531
1532 notes.push(
1533 "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."
1534 .to_string(),
1535 );
1536
1537 Ok(ConformerSearchResult {
1538 generated: generated_count,
1539 unique: unique.len(),
1540 rotatable_bonds,
1541 conformers: unique,
1542 clusters,
1543 notes,
1544 })
1545}
1546
1547pub fn compute_uff_energy_with_aromatic_heuristics(
1549 smiles: &str,
1550 coords: &[f64],
1551) -> Result<reactivity::UffHeuristicEnergy, String> {
1552 let mol = graph::Molecule::from_smiles(smiles)?;
1553 let n = mol.graph.node_count();
1554 if coords.len() != n * 3 {
1555 return Err(format!("coords length {} != 3 * atoms {}", coords.len(), n));
1556 }
1557
1558 let ff = forcefield::builder::build_uff_force_field(&mol);
1559 let mut gradient = vec![0.0f64; n * 3];
1560 let raw = ff.compute_system_energy_and_gradients(coords, &mut gradient);
1561 Ok(reactivity::apply_aromatic_uff_correction(&mol, raw))
1562}
1563
1564pub fn compute_empirical_pka(smiles: &str) -> Result<reactivity::EmpiricalPkaResult, String> {
1566 let mol = graph::Molecule::from_smiles(smiles)?;
1567 let charges = compute_charges(smiles)?;
1568 Ok(reactivity::estimate_empirical_pka(&mol, &charges.charges))
1569}
1570
1571pub fn create_unit_cell(
1573 a: f64,
1574 b: f64,
1575 c: f64,
1576 alpha: f64,
1577 beta: f64,
1578 gamma: f64,
1579) -> materials::UnitCell {
1580 materials::UnitCell::from_parameters(&materials::CellParameters {
1581 a,
1582 b,
1583 c,
1584 alpha,
1585 beta,
1586 gamma,
1587 })
1588}
1589
1590pub fn assemble_framework(
1594 node: &materials::Sbu,
1595 linker: &materials::Sbu,
1596 topology: &materials::Topology,
1597 cell: &materials::UnitCell,
1598) -> materials::CrystalStructure {
1599 materials::assemble_framework(node, linker, topology, cell)
1600}
1601
1602#[cfg(test)]
1603mod tests {
1604 use super::*;
1605
1606 #[test]
1607 fn test_cisplatin_style_support_metadata() {
1608 let smiles = "[Pt](Cl)(Cl)([NH3])[NH3]";
1609 let mol = parse(smiles).expect("Cisplatin-style example should parse");
1610 let elements: Vec<u8> = (0..mol.graph.node_count())
1611 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
1612 .collect();
1613
1614 let caps = get_system_capabilities(&elements);
1615 assert!(caps.eht.available);
1616 assert!(matches!(
1617 caps.eht.confidence,
1618 eht::SupportLevel::Experimental
1619 ));
1620 }
1621
1622 #[test]
1623 fn test_pt_system_routes_to_uff_when_experimental_disabled() {
1624 let smiles = "[Pt](Cl)(Cl)([NH3])[NH3]";
1625 let mol = parse(smiles).expect("Pt example should parse");
1626 let n = mol.graph.node_count();
1627 let elements: Vec<u8> = (0..n)
1628 .map(|i| mol.graph[petgraph::graph::NodeIndex::new(i)].element)
1629 .collect();
1630 let positions = vec![[0.0, 0.0, 0.0]; n];
1631
1632 let result = compute_eht_or_uff_fallback(smiles, &elements, &positions, false)
1633 .expect("Fallback workflow should return a result");
1634
1635 assert!(matches!(
1636 result,
1637 ElectronicWorkflowResult::UffFallback { .. }
1638 ));
1639 }
1640
1641 #[test]
1642 fn test_method_plan_prefers_supported_workflows_for_organic_systems() {
1643 let plan = get_system_method_plan(&[6, 1, 1, 1, 1]);
1644
1645 assert_eq!(plan.geometry.recommended, Some(ScientificMethod::Embed));
1646 assert_eq!(
1647 plan.force_field_energy.recommended,
1648 Some(ScientificMethod::Uff)
1649 );
1650 assert_eq!(plan.orbitals.recommended, Some(ScientificMethod::Eht));
1651 assert_eq!(plan.population.recommended, Some(ScientificMethod::Eht));
1652 assert_eq!(plan.orbital_grid.recommended, Some(ScientificMethod::Eht));
1653 assert!(plan.orbitals.methods[0].confidence_score > 0.9);
1654 }
1655
1656 #[test]
1657 fn test_method_plan_marks_metal_orbital_workflow_experimental() {
1658 let plan = get_system_method_plan(&[78, 17, 17, 7, 7]);
1659
1660 assert_eq!(plan.orbitals.recommended, Some(ScientificMethod::Eht));
1661 assert_eq!(
1662 plan.force_field_energy.recommended,
1663 Some(ScientificMethod::Uff)
1664 );
1665 assert!(matches!(
1666 plan.orbitals.methods[0].confidence,
1667 eht::SupportLevel::Experimental
1668 ));
1669 assert!(plan.orbitals.methods[0].confidence_score < 0.9);
1670 assert!(!plan.orbitals.methods[0].warnings.is_empty());
1671 }
1672
1673 #[test]
1674 fn test_method_plan_reports_unavailable_workflows_for_unsupported_elements() {
1675 let plan = get_system_method_plan(&[92]);
1676
1677 assert_eq!(plan.force_field_energy.recommended, None);
1678 assert_eq!(plan.orbitals.recommended, None);
1679 assert_eq!(plan.population.recommended, None);
1680 assert_eq!(plan.orbital_grid.recommended, None);
1681 assert!(!plan.orbitals.methods[0].limitations.is_empty());
1682 }
1683
1684 #[test]
1685 fn test_compare_methods_supported_system_returns_success_rows() {
1686 let result = compare_methods("CC", &[6, 6], &[[0.0, 0.0, 0.0], [1.54, 0.0, 0.0]], true);
1687 assert_eq!(result.comparisons.len(), 2);
1688 assert!(result
1689 .comparisons
1690 .iter()
1691 .any(|entry| matches!(entry.method, ScientificMethod::Uff) && entry.available));
1692 assert!(result.comparisons.iter().any(|entry| matches!(
1693 entry.method,
1694 ScientificMethod::Eht
1695 ) && matches!(
1696 entry.status,
1697 MethodComparisonStatus::Success
1698 )));
1699 }
1700
1701 #[test]
1702 fn test_compare_methods_blocks_experimental_eht_when_disabled() {
1703 let result = compare_methods("[O]", &[78], &[[0.0, 0.0, 0.0]], false);
1704 let eht_row = result
1705 .comparisons
1706 .iter()
1707 .find(|entry| matches!(entry.method, ScientificMethod::Eht))
1708 .expect("EHT row must exist");
1709 assert!(matches!(
1710 eht_row.status,
1711 MethodComparisonStatus::Unavailable
1712 ));
1713 }
1714
1715 #[test]
1716 fn test_compare_methods_reports_unavailable_for_unsupported_elements() {
1717 let result = compare_methods("[O]", &[92], &[[0.0, 0.0, 0.0]], true);
1718 assert!(result
1719 .comparisons
1720 .iter()
1721 .all(|entry| matches!(entry.status, MethodComparisonStatus::Unavailable)));
1722 }
1723
1724 #[test]
1725 fn test_compute_fukui_descriptors_returns_atomwise_output() {
1726 let elements = [8u8, 1, 1];
1727 let positions = [[0.0, 0.0, 0.0], [0.757, 0.586, 0.0], [-0.757, 0.586, 0.0]];
1728 let result = compute_fukui_descriptors(&elements, &positions).unwrap();
1729 assert_eq!(result.num_atoms, 3);
1730 assert_eq!(result.condensed.len(), 3);
1731 }
1732
1733 #[test]
1734 fn test_compute_uv_vis_spectrum_returns_requested_grid_size() {
1735 let elements = [6u8, 6, 1, 1, 1, 1];
1736 let positions = [
1737 [0.0, 0.0, 0.0],
1738 [1.34, 0.0, 0.0],
1739 [-0.6, 0.92, 0.0],
1740 [-0.6, -0.92, 0.0],
1741 [1.94, 0.92, 0.0],
1742 [1.94, -0.92, 0.0],
1743 ];
1744 let spectrum = compute_uv_vis_spectrum(&elements, &positions, 0.2, 0.5, 8.0, 256).unwrap();
1745 assert_eq!(spectrum.energies_ev.len(), 256);
1746 assert_eq!(spectrum.intensities.len(), 256);
1747 }
1748
1749 #[test]
1750 fn test_analyze_graph_features_reports_benzene_aromaticity() {
1751 let analysis = analyze_graph_features("c1ccccc1").unwrap();
1752 assert_eq!(analysis.aromaticity.aromatic_atoms.len(), 12);
1753 assert_eq!(
1754 analysis
1755 .aromaticity
1756 .aromatic_atoms
1757 .iter()
1758 .filter(|v| **v)
1759 .count(),
1760 6
1761 );
1762 assert_eq!(analysis.aromaticity.aromatic_bonds.len(), 6);
1763 }
1764
1765 #[test]
1766 fn test_compute_empirical_pka_finds_acidic_site_for_acetic_acid() {
1767 let result = compute_empirical_pka("CC(=O)O").unwrap();
1768 assert!(!result.acidic_sites.is_empty());
1769 }
1770
1771 #[test]
1772 fn test_compute_uff_energy_with_aromatic_heuristics_applies_correction() {
1773 let conf = embed("c1ccccc1", 42);
1774 assert!(conf.error.is_none());
1775
1776 let result = compute_uff_energy_with_aromatic_heuristics("c1ccccc1", &conf.coords).unwrap();
1777 assert!(result.aromatic_bond_count >= 6);
1778 assert!(result.corrected_energy_kcal_mol <= result.raw_energy_kcal_mol);
1779 }
1780
1781 #[test]
1782 fn test_search_conformers_with_uff_returns_ranked_unique_ensemble() {
1783 let result = search_conformers_with_uff("CCCC", 10, 42, 0.2).unwrap();
1784 assert!(result.generated >= 1);
1785 assert!(result.unique >= 1);
1786 assert_eq!(result.unique, result.conformers.len());
1787 assert_eq!(result.unique, result.clusters.len());
1788 assert!(result.rotatable_bonds >= 1);
1789
1790 let mut total_members = 0usize;
1791 for (i, cluster) in result.clusters.iter().enumerate() {
1792 assert_eq!(cluster.cluster_id, i);
1793 assert!(cluster.size >= 1);
1794 total_members += cluster.size;
1795 assert_eq!(result.conformers[i].cluster_id, Some(i));
1796 assert_eq!(result.conformers[i].seed, cluster.representative_seed);
1797 }
1798 assert_eq!(total_members, result.generated);
1799
1800 for i in 1..result.conformers.len() {
1801 assert!(
1802 result.conformers[i - 1].energy_kcal_mol <= result.conformers[i].energy_kcal_mol
1803 );
1804 }
1805 }
1806
1807 #[test]
1808 fn test_search_conformers_with_uff_large_rmsd_threshold_collapses_duplicates() {
1809 let result = search_conformers_with_uff("CCCC", 8, 123, 10.0).unwrap();
1810 assert_eq!(result.unique, 1);
1811 assert_eq!(result.clusters.len(), 1);
1812 assert_eq!(result.clusters[0].size, result.generated);
1813 }
1814
1815 #[test]
1816 fn test_compute_md_trajectory_velocity_verlet_runs() {
1817 let conf = embed("CC", 42);
1818 assert!(conf.error.is_none());
1819
1820 let trj = compute_md_trajectory("CC", &conf.coords, 10, 0.25, 7).unwrap();
1821 assert_eq!(trj.frames.len(), 11);
1822 assert!(trj
1823 .frames
1824 .iter()
1825 .all(|f| f.coords.iter().all(|v| v.is_finite())));
1826 }
1827
1828 #[test]
1829 fn test_compute_md_trajectory_nvt_runs() {
1830 let conf = embed("CCO", 42);
1831 assert!(conf.error.is_none());
1832
1833 let trj =
1834 compute_md_trajectory_nvt("CCO", &conf.coords, 12, 0.25, 17, 300.0, 10.0).unwrap();
1835 assert_eq!(trj.frames.len(), 13);
1836 assert!(trj.frames.iter().all(|f| f.temperature_k.is_finite()));
1837 }
1838
1839 #[test]
1840 fn test_compute_simplified_neb_path_runs() {
1841 let c1 = embed("CC", 42);
1842 let c2 = embed("CC", 43);
1843 assert!(c1.error.is_none());
1844 assert!(c2.error.is_none());
1845
1846 let path =
1847 compute_simplified_neb_path("CC", &c1.coords, &c2.coords, 6, 20, 0.01, 1e-5).unwrap();
1848 assert_eq!(path.images.len(), 6);
1849 assert!(path
1850 .images
1851 .iter()
1852 .all(|img| img.potential_energy_kcal_mol.is_finite()));
1853 }
1854}